32 #include <CLHEP/Units/SystemOfUnits.h> 33 #include <CLHEP/Units/GlobalPhysicalConstants.h> 40 for (
int i = 0;
i < 6; ++
i)
68 void fillHitsInfo(std::pair<hitsinfo, energysum> hit_,
unsigned int itimeslice,
double esum);
73 const char* yTitle =
"",
79 const char* yTitle =
"",
108 : nameDetector_(iConfig.getParameter<
std::
string>(
"DetectorName")),
109 caloHitSource_(iConfig.getParameter<
std::
string>(
"CaloHitSource")),
110 times_(iConfig.getParameter<
std::
vector<double> >(
"TimeSlices")),
111 verbosity_(iConfig.getUntrackedParameter<
int>(
"Verbosity", 0)),
114 tok_hepMC_(consumes<edm::HepMCProduct>(
edm::InputTag(
"generatorSmeared"))),
115 tok_hits_(consumes<edm::PCaloHitContainer>(
edm::InputTag(
"g4SimHits", caloHitSource_))),
117 nTimes_ = (times_.size() > maxTime_) ? maxTime_ : times_.size();
122 std::vector<double> times = {25.0, 1000.0};
125 desc.add<std::vector<double> >(
"TimeSlices", times);
126 desc.addUntracked<
int>(
"Verbosity", 0);
127 desc.addUntracked<
bool>(
"TestNumber",
true);
128 descriptions.
add(
"hgcalSimHitValidationEE",
desc);
140 for (HepMC::GenEvent::particle_const_iterator
p = myGenEvent->particles_begin();
p != myGenEvent->particles_end();
142 edm::LogVerbatim(
"HGCalValidation") <<
"Particle[" <<
k <<
"] with pt " << (*p)->momentum().perp() <<
" eta " 143 << (*p)->momentum().eta() <<
" phi " << (*p)->momentum().phi();
150 if (theCaloHitContainers.
isValid()) {
152 edm::LogVerbatim(
"HGCalValidation") <<
" PcalohitItr = " << theCaloHitContainers->size();
153 std::vector<PCaloHit> caloHits;
154 caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
162 std::map<int, int> OccupancyMap_plus, OccupancyMap_minus;
163 OccupancyMap_plus.clear();
164 OccupancyMap_minus.clear();
166 std::map<uint32_t, std::pair<hitsinfo, energysum> > map_hits;
171 unsigned int nused(0);
172 for (
unsigned int i = 0;
i <
hits.size();
i++) {
175 uint32_t id_ =
hits[
i].id();
180 cell =
detId.cellU();
203 <<
" type = " <<
type <<
" layer = " <<
layer <<
" cell = " << cell <<
":" <<
cell2 <<
" energy = " <<
energy 204 <<
" energyem = " <<
hits[
i].energyEM() <<
" energyhad = " <<
hits[
i].energyHad() <<
" time = " <<
time;
206 HepGeom::Point3D<float> gcoord;
207 std::pair<float, float>
xy;
209 xy =
hgcons_->
locateCell(
zside,
layer,
sector,
subsector, cell,
cell2,
false,
true,
false,
false,
false);
216 float xp = (zp < 0) ? -
xy.first :
xy.first;
217 gcoord = HepGeom::Point3D<float>(xp,
xy.second, zp);
218 double tof = (gcoord.mag() * CLHEP::mm) / CLHEP::c_light;
221 << std::hex << id_ <<
std::dec <<
" global coordinate " << gcoord <<
" time " <<
time <<
":" << tof;
226 if (map_hits.count(id_) != 0) {
227 hinfo = map_hits[id_].first;
228 esum = map_hits[id_].second;
230 hinfo.x = gcoord.x();
231 hinfo.y = gcoord.y();
232 hinfo.z = gcoord.z();
239 hinfo.phi = gcoord.getPhi();
240 hinfo.eta = gcoord.getEta();
250 <<
" gz = " <<
hinfo.z <<
" phi = " <<
hinfo.phi <<
" eta = " <<
hinfo.eta;
251 map_hits[id_] = std::pair<hitsinfo, energysum>(
hinfo, esum);
255 <<
" detector elements being hit";
257 std::map<uint32_t, std::pair<hitsinfo, energysum> >::iterator itr;
258 for (itr = map_hits.begin(); itr != map_hits.end(); ++itr) {
264 int partialType = -1;
271 for (
unsigned int itimeslice = 0; itimeslice <
nTimes_; itimeslice++) {
283 edm::LogVerbatim(
"HGCalValidation") <<
"With map:used:total " <<
hits.size() <<
"|" << nused <<
"|" 284 << map_hits.size() <<
" hits";
286 for (
auto const& itr : OccupancyMap_plus) {
287 int layer = itr.first;
288 int occupancy = itr.second;
291 for (
auto const& itr : OccupancyMap_minus) {
292 int layer = itr.first;
293 int occupancy = itr.second;
299 if (OccupancyMap.find(
layer) != OccupancyMap.end()) {
300 ++OccupancyMap[
layer];
302 OccupancyMap[
layer] = 1;
307 unsigned int ilayer =
hits.first.layer;
309 energy_[itimeslice].at(ilayer)->Fill(esum);
310 if (itimeslice == 0) {
317 <<
"Problematic Hit for " <<
nameDetector_ <<
" at sector " <<
hits.first.sector <<
":" <<
hits.first.sector2
318 <<
" layer " <<
hits.first.layer <<
" cell " <<
hits.first.cell <<
":" <<
hits.first.cell2 <<
" energy " 319 <<
hits.second.etotal;
327 esum.
eTime[0] * CLHEP::GeV /
330 unsigned int ilayer =
hinfo.layer;
331 double x =
hinfo.x * CLHEP::mm / CLHEP::cm;
332 double y =
hinfo.y * CLHEP::mm / CLHEP::cm;
336 if (!TMath::AreEqualAbs(edep, 0.0, 1.
e-5)) {
338 if (partialType == 0)
344 if (partialType == 0)
350 if (partialType == 0)
387 std::ostringstream histoname;
388 for (
unsigned int il = 0; il <
layers_; ++il) {
391 while (istr1.size() < 2) {
392 istr1.insert(0,
"0");
395 histoname <<
"HitOccupancy_Plus_layer_" << istr1;
398 histoname <<
"HitOccupancy_Minus_layer_" << istr1;
402 histoname <<
"EtaPhi_Plus_" 403 <<
"layer_" << istr1;
406 histoname <<
"EtaPhi_Minus_" 407 <<
"layer_" << istr1;
411 for (
unsigned int itimeslice = 0; itimeslice <
nTimes_; itimeslice++) {
413 histoname <<
"energy_time_" << itimeslice <<
"_layer_" << istr1;
414 energy_[itimeslice].push_back(iB.
book1D(histoname.str().c_str(),
"energy_", 100, 0, 0.1));
420 histoname <<
"energy_FullWafer_Fine_layer_" << istr1;
421 TH1F* hEdepFWF =
createHisto(histoname.str(), 100, 0., 400.,
false);
422 histoSetting(hEdepFWF,
"Eloss (keV)",
"", kRed, kRed, 2);
427 histoname <<
"energy_FullWafer_CoarseThin_layer_" << istr1;
428 TH1F* hEdepFWCN =
createHisto(histoname.str(), 100, 0., 400.,
false);
429 histoSetting(hEdepFWCN,
"Eloss (keV)",
"", kGreen + 1, kGreen + 1, 2);
434 histoname <<
"energy_FullWafer_CoarseThick_layer_" << istr1;
435 TH1F* hEdepFWCK =
createHisto(histoname.str(), 100, 0., 400.,
false);
436 histoSetting(hEdepFWCK,
"Eloss (keV)",
"", kMagenta, kMagenta, 2);
446 histoname <<
"energy_PartialWafer_Fine_layer_" << istr1;
447 TH1F* hEdepPWF =
createHisto(histoname.str(), 100, 0., 400.,
false);
448 histoSetting(hEdepPWF,
"Eloss (keV)",
"", kRed, kRed, 2);
453 histoname <<
"energy_PartialWafer_CoarseThin_layer_" << istr1;
454 TH1F* hEdepPWCN =
createHisto(histoname.str(), 100, 0., 400.,
false);
455 histoSetting(hEdepPWCN,
"Eloss (keV)",
"", kGreen + 1, kGreen + 1, 2);
460 histoname <<
"energy_PartialWafer_CoarseThick_layer_" << istr1;
461 TH1F* hEdepPWCK =
createHisto(histoname.str(), 100, 0., 400.,
false);
462 histoSetting(hEdepPWCK,
"Eloss (keV)",
"", kMagenta, kMagenta, 2);
471 histoname <<
"hitXY_FullWafer_Fine_layer_" << istr1;
472 TH2F* hitXYFWF =
new TH2F(
473 Form(
"hitXYFWF_%s", histoname.str().c_str()), histoname.str().c_str(), 100, -300., 300., 100, -300., 300.);
479 histoname <<
"hitXY_FullWafer_CoarseThin_layer_" << istr1;
480 TH2F* hitXYFWCN =
new TH2F(
481 Form(
"hitXYFWCN_%s", histoname.str().c_str()), histoname.str().c_str(), 100, -300., 300., 100, -300., 300.);
482 histoSetting(hitXYFWCN,
"x (cm)",
"y (cm)", kGreen + 1, kGreen + 1);
487 histoname <<
"hitXY_FullWafer_CoarseThick_layer_" << istr1;
488 TH2F* hitXYFWCK =
new TH2F(
489 Form(
"hitXYFWCK_%s", histoname.str().c_str()), histoname.str().c_str(), 100, -300., 300., 100, -300., 300.);
490 histoSetting(hitXYFWCK,
"x (cm)",
"y (cm)", kMagenta, kMagenta);
497 histoname <<
"hitXY_Scintillator_layer_" << istr1;
498 TH2F* hitXYB =
new TH2F(
499 Form(
"hitXYB_%s", histoname.str().c_str()), histoname.str().c_str(), 100, -300., 300., 100, -300., 300.);
501 hitXYB_.push_back(iB.
book2D(histoname.str().c_str(), hitXYB));
511 std::string histname,
const int nbins,
float minIndexX,
float maxIndexX,
bool isLogX) {
512 TH1F*
hist =
nullptr;
515 double dx = (maxIndexX - minIndexX) /
nbins;
517 xbins[
i] = TMath::Power(10, (minIndexX +
i *
dx));
519 hist =
new TH1F(Form(
"hEdep_%s", histname.c_str()), histname.c_str(),
nbins,
xbins);
521 hist =
new TH1F(Form(
"hEdep_%s", histname.c_str()), histname.c_str(),
nbins, minIndexX, maxIndexX);
532 histo->GetXaxis()->SetTitle(xTitle);
533 histo->GetYaxis()->SetTitle(yTitle);
542 histo->SetMarkerStyle(kFullCircle);
543 histo->SetMarkerSize(0.6);
544 histo->GetXaxis()->SetTitle(xTitle);
545 histo->GetYaxis()->SetTitle(yTitle);
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * MeanHitOccupancy_Plus_
double waferZ(int layer, bool reco) const
Log< level::Info, true > LogVerbatim
void fillOccupancyMap(std::map< int, int > &OccupancyMap, int layer)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
std::vector< MonitorElement * > HitOccupancy_Plus_
void analyze(const edm::Event &, const edm::EventSetup &) override
const std::vector< double > times_
virtual void setCurrentFolder(std::string const &fullpath)
bool waferHexagon8() const
void analyzeHits(std::vector< PCaloHit > &hits)
std::vector< MonitorElement * > energyPWF_
std::vector< MonitorElement * > HitOccupancy_Minus_
Log< level::Error, false > LogError
HGCalGeometryMode::GeometryMode geomMode() const
static std::string to_string(const XMLCh *ch)
MonitorElement * MeanHitOccupancy_Minus_
std::vector< MonitorElement * > energyPWCN_
const std::string nameDetector_
std::vector< MonitorElement * > hitXYFWCK_
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
static const unsigned int maxTime_
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
TH1F * createHisto(std::string histname, const int nbins, float minIndexX, float maxIndexX, bool isLogX=true)
#define DEFINE_FWK_MODULE(type)
unsigned int layers(bool reco) const
std::vector< MonitorElement * > energyFWF_
bool tileTrapezoid() const
std::vector< MonitorElement * > EtaPhi_Plus_
std::pair< float, float > locateCellTrap(int zside, int lay, int ieta, int iphi, bool reco, bool debug) const
std::vector< MonitorElement * > energyFWCN_
const HepMC::GenEvent * GetEvent() const
HGCalSimHitValidation(const edm::ParameterSet &)
~HGCalSimHitValidation() override=default
const edm::EDGetTokenT< edm::HepMCProduct > tok_hepMC_
std::vector< MonitorElement * > energyFWCK_
std::vector< MonitorElement * > hitXYFWF_
const std::string caloHitSource_
std::vector< MonitorElement * > energyPWCK_
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
int waferType(DetId const &id, bool fromFile) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
std::vector< MonitorElement * > energy_[maxTime_]
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void histoSetting(TH1F *&histo, const char *xTitle, const char *yTitle="", Color_t lineColor=kBlack, Color_t markerColor=kBlack, int linewidth=1)
const HGCalDDDConstants * hgcons_
std::vector< MonitorElement * > EtaPhi_Minus_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
std::vector< MonitorElement * > hitXYFWCN_
void fillMuonTomoHistos(int partialType, std::pair< hitsinfo, energysum > hit_)
std::vector< MonitorElement * > hitXYB_
void fillHitsInfo(std::pair< hitsinfo, energysum > hit_, unsigned int itimeslice, double esum)
const edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > tok_hgcal_