CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SimG4HcalValidation.cc
Go to the documentation of this file.
1 // File: SimG4Validation.cc
3 // Description: Main analysis class for Hcal Validation of G4 Hits
11 
12 // to retreive hits
18 
23 
28 
32 
35 
36 #include "G4HCofThisEvent.hh"
37 #include "G4SDManager.hh"
38 #include "G4Step.hh"
39 #include "G4Track.hh"
40 #include "G4VProcess.hh"
41 
42 #include <iostream>
43 #include <memory>
44 #include <string>
45 #include <vector>
46 
47 using namespace geant_units::operators;
48 
50  public Observer<const BeginOfRun *>,
51  public Observer<const BeginOfEvent *>,
52  public Observer<const EndOfEvent *>,
53  public Observer<const G4Step *> {
54 public:
56  SimG4HcalValidation(const SimG4HcalValidation &) = delete; // stop default
57  const SimG4HcalValidation &operator=(const SimG4HcalValidation &) = delete;
58  ~SimG4HcalValidation() override;
59 
60  void registerConsumes(edm::ConsumesCollector) override;
61  void produce(edm::Event &, const edm::EventSetup &) override;
62  void beginRun(edm::EventSetup const &) override;
63 
64 private:
65  void init();
66 
67  // observer classes
68  void update(const BeginOfRun *run) override;
69  void update(const BeginOfEvent *evt) override;
70  void update(const G4Step *step) override;
71  void update(const EndOfEvent *evt) override;
72 
73  // jetfinding and analysis-related stuff
74  void fill(const EndOfEvent *ev);
75  void layerAnalysis(PHcalValidInfoLayer &);
76  void nxNAnalysis(PHcalValidInfoNxN &);
77  void jetAnalysis(PHcalValidInfoJets &);
78  void fetchHits(PHcalValidInfoLayer &);
79  void clear();
80  void collectEnergyRdir(const double, const double);
81  double getHcalScale(std::string, int) const;
82 
83 private:
85 
86  // Keep parameters to instantiate Jet finder later
87  std::unique_ptr<SimG4HcalHitJetFinder> jetf;
88 
89  // Keep reference to instantiate HcalNumberingFromDDD later
90  std::unique_ptr<HcalNumberingFromDDD> numberingFromDDD;
91 
92  // Keep parameters to instantiate HcalTestNumberingScheme later
93  std::unique_ptr<HcalTestNumberingScheme> org;
94 
95  // Hit cache for cluster analysis
96  std::vector<CaloHit> hitcache; // e, eta, phi, time, layer, calo type
97 
98  // scale factors :
99  std::vector<float> scaleHB;
100  std::vector<float> scaleHE;
101  std::vector<float> scaleHF;
102 
103  // to read from parameter set
104  std::vector<std::string> names;
105  double coneSize, ehitThreshold, hhitThreshold;
106  float timeLowlim, timeUplim, eta0, phi0, jetThreshold;
107  bool applySampling, hcalOnly;
109  std::string labelLayer, labelNxN, labelJets;
110 
111  // eta and phi size of windows around eta0, phi0
112  std::vector<double> dEta;
113  std::vector<double> dPhi;
114 
115  // some private members for ananlysis
116  unsigned int count;
117  double edepEB, edepEE, edepHB, edepHE, edepHO;
118  double edepd[5], edepl[20];
119  double een, hen, hoen; // Energy sum in ECAL, HCAL, HO
120  double vhitec, vhithc, enEcal, enHcal;
121 };
122 
124  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("SimG4HcalValidation");
125  infolevel = m_Anal.getParameter<int>("InfoLevel");
126  hcalOnly = m_Anal.getParameter<bool>("HcalClusterOnly");
127  applySampling = m_Anal.getParameter<bool>("HcalSampling");
128  coneSize = m_Anal.getParameter<double>("ConeSize");
129  ehitThreshold = m_Anal.getParameter<double>("EcalHitThreshold");
130  hhitThreshold = m_Anal.getParameter<double>("HcalHitThreshold");
131  timeLowlim = m_Anal.getParameter<double>("TimeLowLimit");
132  timeUplim = m_Anal.getParameter<double>("TimeUpLimit");
133  jetThreshold = m_Anal.getParameter<double>("JetThreshold");
134  eta0 = m_Anal.getParameter<double>("Eta0");
135  phi0 = m_Anal.getParameter<double>("Phi0");
136  names = m_Anal.getParameter<std::vector<std::string>>("Names");
137  labelLayer = m_Anal.getUntrackedParameter<std::string>("LabelLayerInfo", "HcalInfoLayer");
138  labelNxN = m_Anal.getUntrackedParameter<std::string>("LabelNxNInfo", "HcalInfoNxN");
139  labelJets = m_Anal.getUntrackedParameter<std::string>("LabelJetsInfo", "HcalInfoJets");
140 
141  produces<PHcalValidInfoLayer>(labelLayer);
142  if (infolevel > 0)
143  produces<PHcalValidInfoNxN>(labelNxN);
144  if (infolevel > 1)
145  produces<PHcalValidInfoJets>(labelJets);
146 
147  edm::LogVerbatim("ValidHcal") << "HcalTestAnalysis:: Initialized as observer of begin/end events and "
148  << "of G4step with Parameter values: \n\tInfoLevel = " << infolevel
149  << "\n\thcalOnly = " << hcalOnly << "\n\tapplySampling = " << applySampling
150  << "\n\tconeSize = " << coneSize << "\n\tehitThreshold = " << ehitThreshold
151  << "\n\thhitThreshold = " << hhitThreshold << "\n\tttimeLowlim = " << timeLowlim
152  << "\n\tttimeUplim = " << timeUplim << "\n\teta0 = " << eta0
153  << "\n\tphi0 = " << phi0 << "\nLabels (Layer): " << labelLayer
154  << " (NxN): " << labelNxN << " (Jets): " << labelJets;
155 
156  init();
157 }
158 
160  edm::LogVerbatim("ValidHcal") << "\n --------> Total number of selected entries"
161  << " : " << count;
162 }
163 
166  edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::Initialize ESGetToken for HcalDDDSimConstants";
167 }
168 
170  std::unique_ptr<PHcalValidInfoLayer> productLayer(new PHcalValidInfoLayer);
171  layerAnalysis(*productLayer);
172  e.put(std::move(productLayer), labelLayer);
173 
174  if (infolevel > 0) {
175  std::unique_ptr<PHcalValidInfoNxN> productNxN(new PHcalValidInfoNxN);
176  nxNAnalysis(*productNxN);
177  e.put(std::move(productNxN), labelNxN);
178  }
179 
180  if (infolevel > 1) {
181  std::unique_ptr<PHcalValidInfoJets> productJets(new PHcalValidInfoJets);
182  jetAnalysis(*productJets);
183  e.put(std::move(productJets), labelJets);
184  }
185 }
186 
187 //==================================================================== per RUN
189  float sHB[4] = {117., 117., 178., 217.};
190  float sHE[3] = {178., 178., 178.};
191  float sHF[3] = {2.84, 2.09, 0.}; //
192 
193  float deta[4] = {0.0435, 0.1305, 0.2175, 0.3045};
194  float dphi[4] = {0.0436, 0.1309, 0.2182, 0.3054};
195 
196  int i = 0;
197  for (i = 0; i < 4; i++) {
198  scaleHB.push_back(sHB[i]);
199  }
200  for (i = 0; i < 3; i++) {
201  scaleHE.push_back(sHE[i]);
202  }
203  for (i = 0; i < 3; i++) {
204  scaleHF.push_back(sHF[i]);
205  }
206 
207  // window steps;
208  for (i = 0; i < 4; i++) {
209  dEta.push_back(deta[i]);
210  dPhi.push_back(dphi[i]);
211  }
212 
213  // jetfinder conse size setting
214  jetf = std::make_unique<SimG4HcalHitJetFinder>(coneSize);
215 
216  // counter
217  count = 0;
218 }
219 
221  // Numbering From DDD
222  const HcalDDDSimConstants *hcons = &es.getData(ddconsToken_);
223  edm::LogVerbatim("ValidHcal") << "HcalTestAnalysis:: Initialise HcalNumberingFromDDD";
224  numberingFromDDD = std::make_unique<HcalNumberingFromDDD>(hcons);
225 
226  // Numbering scheme
227  org = std::make_unique<HcalTestNumberingScheme>(false);
228 }
229 
231  int irun = (*run)()->GetRunID();
232 
233  edm::LogVerbatim("ValidHcal") << " =====> Begin of Run = " << irun;
234 
235  std::string sdname = names[0];
236  G4SDManager *sd = G4SDManager::GetSDMpointerIfExist();
237  if (sd != nullptr) {
238  G4VSensitiveDetector *aSD = sd->FindSensitiveDetector(sdname);
239  if (aSD == nullptr) {
240  edm::LogWarning("ValidHcal") << "SimG4HcalValidation::beginOfRun: No SD"
241  << " with name " << sdname << " in this "
242  << "Setup";
243  } else {
244  HCalSD *theCaloSD = dynamic_cast<HCalSD *>(aSD);
245  edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::beginOfRun: Finds SD with name " << theCaloSD->GetName()
246  << " in this Setup";
247  if (org.get()) {
248  theCaloSD->setNumberingScheme(org.get());
249  edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::beginOfRun: set a new numbering scheme";
250  }
251  }
252  } else {
253  edm::LogWarning("ValidHcal") << "SimG4HcalValidation::beginOfRun: Could "
254  << "not get SD Manager!";
255  }
256 }
257 
258 //=================================================================== per EVENT
260  int i = 0;
261  edepEB = edepEE = edepHB = edepHE = edepHO = 0.;
262  for (i = 0; i < 5; i++)
263  edepd[i] = 0.;
264  for (i = 0; i < 20; i++)
265  edepl[i] = 0.;
266  vhitec = vhithc = enEcal = enHcal = 0;
267  // Cache reset
268  clear();
269 
270  int iev = (*evt)()->GetEventID();
271  LogDebug("ValidHcal") << "SimG4HcalValidation: =====> Begin of event = " << iev;
272 }
273 
274 //=================================================================== each STEP
275 void SimG4HcalValidation::update(const G4Step *aStep) {
276  if (aStep != nullptr) {
277  G4VPhysicalVolume *curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
278  G4String name = curPV->GetName();
279  name.assign(name, 0, 3);
280  double edeposit = aStep->GetTotalEnergyDeposit();
281  int layer = -1, depth = -1;
282  if (name == "EBR") {
283  depth = 0;
284  edepEB += edeposit;
285  } else if (name == "EFR") {
286  depth = 0;
287  edepEE += edeposit;
288  } else if (name == "HBS") {
289  layer = (curPV->GetCopyNo() / 10) % 100;
290  depth = (curPV->GetCopyNo()) % 10 + 1;
291  if (depth > 0 && depth < 4 && layer >= 0 && layer < 17) {
292  edepHB += edeposit;
293  } else {
294  edm::LogWarning("ValidHcal") << "SimG4HcalValidation:Error " << curPV->GetName() << curPV->GetCopyNo();
295  depth = -1;
296  layer = -1;
297  }
298  } else if (name == "HES") {
299  layer = (curPV->GetCopyNo() / 10) % 100;
300  depth = (curPV->GetCopyNo()) % 10 + 1;
301  if (depth > 0 && depth < 3 && layer >= 0 && layer < 19) {
302  edepHE += edeposit;
303  } else {
304  edm::LogWarning("ValidHcal") << "SimG4HcalValidation:Error " << curPV->GetName() << curPV->GetCopyNo();
305  depth = -1;
306  layer = -1;
307  }
308  } else if (name == "HTS") {
309  layer = (curPV->GetCopyNo() / 10) % 100;
310  depth = (curPV->GetCopyNo()) % 10 + 1;
311  if (depth > 3 && depth < 5 && layer >= 17 && layer < 20) {
312  edepHO += edeposit;
313  } else {
314  edm::LogWarning("ValidHcal") << "SimG4HcalValidation:Error " << curPV->GetName() << curPV->GetCopyNo();
315  depth = -1;
316  layer = -1;
317  }
318  }
319  if (depth >= 0 && depth < 5)
320  edepd[depth] += edeposit;
321  if (layer >= 0 && layer < 20)
322  edepl[layer] += edeposit;
323 
324  if (layer >= 0 && layer < 20) {
325  LogDebug("ValidHcal") << "SimG4HcalValidation:: G4Step: " << name << " Layer " << std::setw(3) << layer
326  << " Depth " << std::setw(2) << depth << " Edep " << std::setw(6) << edeposit / MeV
327  << " MeV";
328  }
329  }
330 }
331 
332 //================================================================ End of EVENT
334  count++;
335 
336  // Fill hits cache for jetfinding etc.
337  fill(evt);
338  LogDebug("ValidHcal") << "SimG4HcalValidation:: --- after Fill ";
339 }
340 
341 //---------------------------------------------------
343  LogDebug("ValidHcal") << "SimG4HcalValidation:Fill event " << (*evt)()->GetEventID();
344 
345  // access to the G4 hit collections
346  G4HCofThisEvent *allHC = (*evt)()->GetHCofThisEvent();
347 
348  int nhc = 0, j = 0;
349 
350  // Hcal
351  int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[0]);
352  CaloG4HitCollection *theHCHC = (CaloG4HitCollection *)allHC->GetHC(HCHCid);
353  LogDebug("ValidHcal") << "SimG4HcalValidation :: Hit Collection for " << names[0] << " of ID " << HCHCid
354  << " is obtained at " << theHCHC;
355  if (HCHCid >= 0 && theHCHC != nullptr) {
356  int nhits = theHCHC->entries();
357  for (j = 0; j < nhits; j++) {
358  CaloG4Hit *aHit = (*theHCHC)[j];
359 
360  double e = aHit->getEnergyDeposit() / GeV;
361  double time = aHit->getTimeSlice();
362 
363  math::XYZPoint pos = aHit->getPosition();
364  double theta = pos.theta();
365  double eta = -log(tan(theta * 0.5));
366  double phi = pos.phi();
367 
368  uint32_t unitID = aHit->getUnitID();
369  int subdet, zside, layer, etaIndex, phiIndex, lay;
370  org->unpackHcalIndex(unitID, subdet, zside, layer, etaIndex, phiIndex, lay);
371 
372  // some logic to separate HO ...
373  layer--;
374  std::string det = "HB";
375  if (subdet == static_cast<int>(HcalForward)) {
376  det = "HF";
377  uint16_t depth = aHit->getDepth();
378  if (depth != 0) {
379  layer += 2;
380  lay += 2;
381  }
382  if (layer == 1)
383  vhithc += e;
384  else if (layer == 0)
385  vhitec += e;
386  else
387  edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::HitPMT " << subdet << " " << (2 * zside - 1) * etaIndex
388  << " " << phiIndex << " " << layer + 1 << " R " << pos.rho() << " Phi "
389  << convertRadToDeg(phi) << " Edep " << e << " Time " << time;
390  } else if (subdet == static_cast<int>(HcalEndcap)) {
391  if (etaIndex <= 20) {
392  det = "HES";
393  } else {
394  det = "HED";
395  }
396  }
397  LogDebug("ValidHcal") << "SimG4HcalValidation::debugFill Hcal " << det << " layer " << std::setw(2) << layer
398  << " lay " << std::setw(2) << lay << " time " << std::setw(6) << time << " theta "
399  << std::setw(8) << theta << " eta " << std::setw(8) << eta << " phi " << std::setw(8) << phi
400  << " e " << std::setw(8) << e << " ID 0x" << std::hex << unitID << " ID dec " << std::dec
401  << (int)unitID;
402 
403  // if desired, apply sampling factors in HCAL !!!
404  if (applySampling)
405  e *= getHcalScale(det, layer);
406 
407  // filter on time & energy
408  if (time >= timeLowlim && time <= timeUplim && e > hhitThreshold) {
409  enHcal += e;
410  CaloHit ahit(subdet, lay, e, eta, phi, time, unitID);
411  hitcache.push_back(ahit); // fill cache
412  ++nhc;
413  }
414  }
415  }
416  LogDebug("ValidHcal") << "SimG4HcalValidation:: HCAL hits : " << nhc;
417 
418  if (!hcalOnly) { //-------------------------- ECAL hits --------------------
419  int ndets = names.size();
420  if (ndets > 3)
421  ndets = 3;
422  for (int idty = 1; idty < ndets; idty++) {
423  std::string det = "EB";
424  if (idty == 2)
425  det = "EE";
426  else if (idty > 2)
427  det = "ES";
428 
429  int nec = 0;
430  int ECHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[idty]);
431  CaloG4HitCollection *theECHC = (CaloG4HitCollection *)allHC->GetHC(ECHCid);
432  LogDebug("ValidHcal") << "SimG4HcalValidation:: Hit Collection for " << names[idty] << " of ID " << ECHCid
433  << " is obtained at " << theECHC;
434  int theechc_entries = theECHC->entries();
435  if (ECHCid >= 0 && theECHC != nullptr) {
436  for (j = 0; j < theechc_entries; j++) {
437  CaloG4Hit *aHit = (*theECHC)[j];
438 
439  double e = aHit->getEnergyDeposit() / GeV;
440  double time = aHit->getTimeSlice();
441 
442  math::XYZPoint pos = aHit->getPosition();
443  double theta = pos.theta();
444  double eta = -log(tan(theta / 2.));
445  double phi = pos.phi();
446  HcalNumberingFromDDD::HcalID id = numberingFromDDD->unitID(eta, phi, 1, 1);
447  uint32_t unitID = org->getUnitID(id);
448  int subdet, zside, layer, ieta, iphi, lay;
449  org->unpackHcalIndex(unitID, subdet, zside, layer, ieta, iphi, lay);
450  subdet = idty + 9;
451  layer = 0;
452  unitID = org->packHcalIndex(subdet, zside, layer, ieta, iphi, lay);
453 
454  // filter on time & energy
455  if (time >= timeLowlim && time <= timeUplim && e > ehitThreshold) {
456  enEcal += e;
457  CaloHit ahit(subdet, lay, e, eta, phi, time, unitID);
458  hitcache.push_back(ahit); // fill cache
459  ++nec;
460  }
461 
462  LogDebug("ValidHcal") << "SimG4HcalValidation::debugFill Ecal " << det << " layer " << std::setw(2) << layer
463  << " lay " << std::setw(2) << lay << " time " << std::setw(6) << time << " theta "
464  << std::setw(8) << theta << " eta " << std::setw(8) << eta << " phi " << std::setw(8)
465  << phi << " e " << std::setw(8) << e << " ID 0x" << std::hex << unitID << " ID dec "
466  << std::dec << (int)unitID;
467  }
468  }
469 
470  LogDebug("ValidHcal") << "SimG4HcalValidation:: " << det << " hits : " << nec;
471  }
472  } // end of if(!hcalOnly)
473 }
474 
476  int i = 0;
477  LogDebug("ValidHcal") << "\n ===>>> SimG4HcalValidation: Energy deposit "
478  << "in MeV\n at EB : " << std::setw(6) << edepEB / MeV << "\n at EE : " << std::setw(6)
479  << edepEE / MeV << "\n at HB : " << std::setw(6) << edepHB / MeV
480  << "\n at HE : " << std::setw(6) << edepHE / MeV << "\n at HO : " << std::setw(6)
481  << edepHO / MeV << "\n ---- SimG4HcalValidation: Energy deposit in";
482  for (i = 0; i < 5; i++)
483  LogDebug("ValidHcal") << " Depth " << std::setw(2) << i << " E " << std::setw(8) << edepd[i] / MeV << " MeV";
484  LogDebug("ValidHcal") << " ---- SimG4HcalValidation: Energy deposit in"
485  << "layers";
486  for (i = 0; i < 20; i++)
487  LogDebug("ValidHcal") << " Layer " << std::setw(2) << i << " E " << std::setw(8) << edepl[i] / MeV << " MeV";
488 
489  product.fillLayers(edepl, edepd, edepHO, edepHB + edepHE, edepEB + edepEE);
490 
491  // Hits in HF
492  product.fillHF(vhitec, vhithc, enEcal, enHcal);
493  LogDebug("ValidHcal") << "SimG4HcalValidation::HF hits " << vhitec << " in EC and " << vhithc << " in HC\n"
494  << " HB/HE " << enEcal << " in EC and " << enHcal << " in HC";
495 
496  // Another HCAL hist to porcess and store separately (a bit more complicated)
497  fetchHits(product);
498 
499  LogDebug("ValidHcal") << " layerAnalysis ===> after fetchHits";
500 }
501 
502 //-----------------------------------------------------------------------------
504  std::vector<CaloHit> *hits = &hitcache;
505  std::vector<CaloHit>::iterator hit_itr;
506 
507  LogDebug("ValidHcal") << "SimG4HcalValidation::NxNAnalysis : entrance ";
508 
509  collectEnergyRdir(eta0, phi0); // HCAL and ECAL energy in SimHitCache
510  // around (eta0,phi0)
511 
512  LogDebug("ValidHcal") << " NxNAnalysis : coolectEnergyRdir - Ecal " << een << " Hcal " << hen;
513 
514  double etot = 0.; // total e deposited in "cluster"
515  double ee = 0.; // ECAL e deposited in "cluster"
516  double he = 0.; // HCAL e deposited in "cluster"
517  double hoe = 0.; // HO e deposited in "cluster"
518 
519  int max = dEta.size(); // 4
520 
521  for (hit_itr = hits->begin(); hit_itr < hits->end(); hit_itr++) {
522  double e = hit_itr->e();
523  double t = hit_itr->t();
524  double eta = hit_itr->eta();
525  double phi = hit_itr->phi();
526  int type = hit_itr->det();
527  int layer = hit_itr->layer();
528 
529  // NxN calulation
530 
531  if (fabs(eta0 - eta) <= dEta[max - 1] && fabs(phi0 - phi) <= dPhi[max - 1]) {
532  etot += e;
533  if (type == 10 || type == 11 || type == 12)
534  ee += e;
535  if (type == static_cast<int>(HcalBarrel) || type == static_cast<int>(HcalEndcap) ||
536  type == static_cast<int>(HcalForward)) {
537  he += e;
538  if (type == static_cast<int>(HcalBarrel) && layer > 17)
539  hoe += e;
540 
541  // which concrete i-th square ?
542  for (int i = 0; i < max; i++) {
543  if ((eta0 - eta) <= dEta[i] && (eta0 - eta) >= -dEta[i] && (phi0 - phi) <= dPhi[i] &&
544  (phi0 - phi) >= -dPhi[i]) {
545  LogDebug("ValidHcal") << "SimG4HcalValidation:: nxNAnalysis eta0,"
546  << " phi0 = " << eta0 << " " << phi0 << " type, layer = " << type << "," << layer
547  << " eta, phi = " << eta << " " << phi;
548 
549  product.fillTProfileNxN(e, i, t);
550  break;
551  }
552  }
553  // here comes break ...
554  }
555  }
556  }
557 
558  product.fillEcollectNxN(ee, he, hoe, etot);
559  product.fillHvsE(een, hen, hoen, een + hen);
560 
561  LogDebug("ValidHcal") << " nxNAnalysis ===> after fillHvsE";
562 }
563 
564 //-----------------------------------------------------------------------------
566  std::vector<CaloHit> *hhit = &hitcache;
567 
568  jetf->setInput(hhit);
569 
570  // zeroing cluster list, perfom clustering, fill cluster list & return pntr
571  std::vector<SimG4HcalHitCluster> *result = jetf->getClusters(hcalOnly);
572 
573  std::vector<SimG4HcalHitCluster>::iterator clus_itr;
574 
575  LogDebug("ValidHcal") << "\n ---------- Final list of " << (*result).size() << " clusters ---------------";
576  for (clus_itr = result->begin(); clus_itr < result->end(); clus_itr++)
577  LogDebug("ValidHcal") << (*clus_itr);
578 
579  std::vector<double> enevec, etavec, phivec;
580 
581  if (!(*result).empty()) {
582  sort((*result).begin(), (*result).end());
583 
584  clus_itr = result->begin(); // first cluster only
585  double etac = clus_itr->eta();
586  double phic = clus_itr->phi();
587 
588  double ecal_collect = 0.; // collect Ecal energy in the cone
589  if (!hcalOnly) {
590  ecal_collect = clus_itr->collectEcalEnergyR();
591  } else {
592  collectEnergyRdir(etac, phic);
593  ecal_collect = een;
594  }
595  LogDebug("ValidHcal") << " JetAnalysis ===> ecal_collect = " << ecal_collect;
596 
597  // eta-phi deviation of the cluster from nominal (eta0, phi0) values
598  double dist = jetf->rDist(eta0, phi0, etac, phic);
599  LogDebug("ValidHcal") << " JetAnalysis ===> eta phi deviation = " << dist;
600  product.fillEtaPhiProfileJet(eta0, phi0, etac, phic, dist);
601 
602  LogDebug("ValidHcal") << " JetAnalysis ===> after fillEtaPhiProfileJet";
603 
604  std::vector<CaloHit> *hits = clus_itr->getHits();
605  std::vector<CaloHit>::iterator hit_itr;
606 
607  double ee = 0., he = 0., hoe = 0., etot = 0.;
608 
609  // cycle over all hits in the FIRST cluster
610  for (hit_itr = hits->begin(); hit_itr < hits->end(); hit_itr++) {
611  double e = hit_itr->e();
612  double t = hit_itr->t();
613  double r = jetf->rDist(&(*clus_itr), &(*hit_itr));
614 
615  // energy collection
616  etot += e;
617  if (hit_itr->det() == 10 || hit_itr->det() == 11 || hit_itr->det() == 12)
618  ee += e;
619  if (hit_itr->det() == static_cast<int>(HcalBarrel) || hit_itr->det() == static_cast<int>(HcalEndcap) ||
620  hit_itr->det() == static_cast<int>(HcalForward)) {
621  he += e;
622  if (hit_itr->det() == static_cast<int>(HcalBarrel) && hit_itr->layer() > 17)
623  hoe += e;
624  }
625 
626  if (hit_itr->det() == static_cast<int>(HcalBarrel) || hit_itr->det() == static_cast<int>(HcalEndcap) ||
627  hit_itr->det() == static_cast<int>(HcalForward)) {
628  product.fillTProfileJet(he, r, t);
629  }
630  }
631 
632  product.fillEcollectJet(ee, he, hoe, etot);
633 
634  LogDebug("ValidHcal") << " JetAnalysis ===> after fillEcollectJet: "
635  << "ee/he/hoe/etot " << ee << "/" << he << "/" << hoe << "/" << etot;
636 
637  // Loop over clusters
638  for (clus_itr = result->begin(); clus_itr < result->end(); clus_itr++) {
639  if (clus_itr->e() > jetThreshold) {
640  enevec.push_back(clus_itr->e());
641  etavec.push_back(clus_itr->eta());
642  phivec.push_back(clus_itr->phi());
643  }
644  }
645  product.fillJets(enevec, etavec, phivec);
646 
647  LogDebug("ValidHcal") << " JetAnalysis ===> after fillJets\n"
648  << " JetAnalysis ===> (*result).size() " << (*result).size();
649 
650  // Di-jet mass
651  if (etavec.size() > 1) {
652  if (etavec[0] > -2.5 && etavec[0] < 2.5 && etavec[1] > -2.5 && etavec[1] < 2.5) {
653  LogDebug("ValidHcal") << " JetAnalysis ===> Di-jet mass enter\n"
654  << " JetAnalysis ===> Di-jet vectors ";
655  for (unsigned int i = 0; i < enevec.size(); i++)
656  LogDebug("ValidHcal") << " e, eta, phi = " << enevec[i] << " " << etavec[i] << " " << phivec[i];
657 
658  double et0 = enevec[0] / cosh(etavec[0]);
659  double px0 = et0 * cos(phivec[0]);
660  double py0 = et0 * sin(phivec[0]);
661  double pz0 = et0 * sinh(etavec[0]);
662  double et1 = enevec[1] / cosh(etavec[1]);
663  double px1 = et1 * cos(phivec[1]);
664  double py1 = et1 * sin(phivec[1]);
665  double pz1 = et1 * sinh(etavec[1]);
666 
667  double dijetmass2 = (enevec[0] + enevec[1]) * (enevec[0] + enevec[1]) - (px1 + px0) * (px1 + px0) -
668  (py1 + py0) * (py1 + py0) - (pz1 + pz0) * (pz1 + pz0);
669 
670  LogDebug("ValidHcal") << " JetAnalysis ===> Di-jet massSQ " << dijetmass2;
671 
672  double dijetmass;
673  if (dijetmass2 >= 0.)
674  dijetmass = sqrt(dijetmass2);
675  else
676  dijetmass = -sqrt(-1. * dijetmass2);
677 
678  product.fillDiJets(dijetmass);
679 
680  LogDebug("ValidHcal") << " JetAnalysis ===> after fillDiJets";
681  }
682  }
683  }
684 }
685 
686 //---------------------------------------------------
688  LogDebug("ValidHcal") << "Enter SimG4HcalValidation::fetchHits with " << hitcache.size() << " hits";
689  int nHit = hitcache.size();
690  int hit = 0;
691  int i;
692  std::vector<CaloHit>::iterator itr;
693  std::vector<CaloHit *> lhits(nHit);
694  for (i = 0, itr = hitcache.begin(); itr != hitcache.end(); i++, itr++) {
695  uint32_t unitID = itr->id();
696  int subdet, zside, group, ieta, iphi, lay;
697  HcalTestNumbering::unpackHcalIndex(unitID, subdet, zside, group, ieta, iphi, lay);
698  subdet = itr->det();
699  lay = itr->layer();
700  group = (subdet & 15) << 20;
701  group += ((lay - 1) & 31) << 15;
702  group += (zside & 1) << 14;
703  group += (ieta & 127) << 7;
704  group += (iphi & 127);
705  itr->setId(group);
706  lhits[i] = &hitcache[i];
707  LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Original " << i << " " << hitcache[i] << "\n"
708  << "SimG4HcalValidation::fetchHits:Copied " << i << " " << *lhits[i];
709  }
710  sort(lhits.begin(), lhits.end(), CaloHitIdMore());
711  std::vector<CaloHit *>::iterator k1, k2;
712  for (i = 0, k1 = lhits.begin(); k1 != lhits.end(); i++, k1++)
713  LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Sorted " << i << " " << **k1;
714  int nHits = 0;
715  for (i = 0, k1 = lhits.begin(); k1 != lhits.end(); i++, k1++) {
716  double ehit = (**k1).e();
717  double t = (**k1).t();
718  uint32_t unitID = (**k1).id();
719  int jump = 0;
720  LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Start " << i << " U/T/E"
721  << " 0x" << std::hex << unitID << std::dec << " " << t << " " << ehit;
722  for (k2 = k1 + 1; k2 != lhits.end() && (t - (**k2).t()) < 1 && (t - (**k2).t()) > -1 && unitID == (**k2).id();
723  k2++) {
724  ehit += (**k2).e();
725  LogDebug("ValidHcal") << "\t + " << (**k2).e();
726  jump++;
727  }
728  LogDebug("ValidHcal") << "\t = " << ehit << " in " << jump;
729 
730  double eta = (*k1)->eta();
731  double phi = (*k1)->phi();
732  int lay = ((unitID >> 15) & 31) + 1;
733  int subdet = (unitID >> 20) & 15;
734  int zside = (unitID >> 14) & 1;
735  int ieta = (unitID >> 7) & 127;
736  int iphi = (unitID)&127;
737 
738  // All hits in cache
739  product.fillHits(nHits, lay, subdet, eta, phi, ehit, t);
740  nHits++;
741 
742  LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Hit " << nHits << " " << i << " ID 0x" << std::hex
743  << unitID << " det " << std::dec << subdet << " " << lay << " " << zside << " " << ieta
744  << " " << iphi << " Time " << t << " E " << ehit;
745 
746  i += jump;
747  k1 += jump;
748  }
749 
750  LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits called with " << nHit << " hits"
751  << " and writes out " << nHits << '(' << hit << ") hits";
752 }
753 //---------------------------------------------------
754 void SimG4HcalValidation::clear() { hitcache.erase(hitcache.begin(), hitcache.end()); }
755 
756 //---------------------------------------------------
757 void SimG4HcalValidation::collectEnergyRdir(const double eta0, const double phi0) {
758  std::vector<CaloHit> *hits = &hitcache;
759  std::vector<CaloHit>::iterator hit_itr;
760 
761  double sume = 0., sumh = 0., sumho = 0.;
762 
763  for (hit_itr = hits->begin(); hit_itr < hits->end(); hit_itr++) {
764  double e = hit_itr->e();
765  double eta = hit_itr->eta();
766  double phi = hit_itr->phi();
767 
768  int type = hit_itr->det();
769 
770  double r = jetf->rDist(eta0, phi0, eta, phi);
771  if (r < coneSize) {
772  if (type == 10 || type == 11 || type == 12) {
773  sume += e;
774  } else {
775  sumh += e;
776  if (type == static_cast<int>(HcalBarrel) && hit_itr->layer() > 17)
777  sumho += e;
778  }
779  }
780  }
781 
782  een = sume;
783  hen = sumh;
784  hoen = sumho;
785 }
786 
787 //---------------------------------------------------
789  double tmp = 0.;
790 
791  if (det == "HB") {
792  tmp = scaleHB[layer];
793  } else if (det == "HES" || det == "HED") {
794  tmp = scaleHE[layer];
795  } else if (det == "HF") {
796  tmp = scaleHF[layer];
797  }
798 
799  return tmp;
800 }
801 
805 
void fillJets(const std::vector< double > &enj, const std::vector< double > &etaj, const std::vector< double > &phij)
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
#define DEFINE_SIMWATCHER(type)
math::XYZPoint getPosition() const
Definition: CaloG4Hit.h:52
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const double GeV
Definition: MathUtil.h:16
static std::vector< std::string > checklist log
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
void fillHvsE(double ee, double he, double hoe, double etot)
void fillTProfileJet(double e, double r, double t)
double getHcalScale(std::string, int) const
void beginRun(edm::EventSetup const &) override
void fill(const EndOfEvent *ev)
void fillLayers(double el[], double ed[], double ho, double hbhe, double ebee)
std::vector< float > scaleHF
void registerConsumes(edm::ConsumesCollector) override
void fillEcollectJet(double ee, double he, double hoe, double etot)
void collectEnergyRdir(const double, const double)
std::vector< std::string > names
void fillEtaPhiProfileJet(double eta0, double phi0, double eta, double phi, double dist)
activeDets clear()
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
std::vector< CaloHit > hitcache
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void produce(edm::Event &, const edm::EventSetup &) override
int init
Definition: HydjetWrapper.h:64
void fillDiJets(double mass)
Geom::Theta< T > theta() const
uint16_t getDepth() const
Definition: CaloG4Hit.h:69
std::vector< double > dEta
bool ev
int zside(DetId const &)
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:544
const std::string names[nVars_]
void fillHits(int Nhits, int lay, int unitID, double eta, double phi, double ehit, double t)
constexpr std::array< uint8_t, layerIndexSize > layer
tuple result
Definition: mps_fire.py:311
bool getData(T &iHolder) const
Definition: EventSetup.h:128
const double MeV
Definition: HCalSD.h:38
void jetAnalysis(PHcalValidInfoJets &)
void nxNAnalysis(PHcalValidInfoNxN &)
std::vector< double > dPhi
void fetchHits(PHcalValidInfoLayer &)
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def move
Definition: eostools.py:511
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
edm::ESGetToken< HcalDDDSimConstants, HcalSimNumberingRecord > ddconsToken_
std::vector< float > scaleHB
void fillEcollectNxN(double een, double hen, double hoen, double etotn)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
tuple group
Definition: watchdog.py:82
void fill(std::map< std::string, TH1 * > &h, const std::string &s, double x)
std::unique_ptr< SimG4HcalHitJetFinder > jetf
void fillTProfileNxN(double e, int i, double t)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double sd
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
std::vector< float > scaleHE
std::unique_ptr< HcalTestNumberingScheme > org
#define update(a, b)
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
double getTimeSlice() const
Definition: CaloG4Hit.h:67
SimG4HcalValidation(const edm::ParameterSet &p)
step
Definition: StallMonitor.cc:94
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
Log< level::Warning, false > LogWarning
HitContainer const *__restrict__ TkSoA const *__restrict__ Quality const *__restrict__ CAHitNtupletGeneratorKernelsGPU::HitToTuple const *__restrict__ int32_t int iev
tmp
align.sh
Definition: createJobs.py:716
void layerAnalysis(PHcalValidInfoLayer &)
void fillHF(double fibl, double fibs, double enec, double enhc)
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
#define LogDebug(id)