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