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