CMS 3D CMS Logo

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