CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HcalTB04Analysis.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestBeam
4 // Class : HcalTB04Analysis
5 //
6 // Implementation:
7 // Main analysis class for Hcal Test Beam 2004 Analysis
8 //
9 // Usage: A Simwatcher class and can be activated from Oscarproducer module
10 //
11 // Original Author:
12 // Created: Tue May 16 10:14:34 CEST 2006
13 //
14 
15 // system include files
16 #include <cmath>
17 #include <iomanip>
18 #include <iostream>
19 #include <memory>
20 #include <vector>
21 #include <string>
22 
23 // user include files
30 
31 // to retreive hits
42 
48 
52 
53 #include "G4SDManager.hh"
54 #include "G4Step.hh"
55 #include "G4Track.hh"
56 #include "G4ThreeVector.hh"
57 #include "G4VProcess.hh"
58 #include "G4HCofThisEvent.hh"
59 
60 #include <CLHEP/Random/RandGaussQ.h>
61 #include <CLHEP/Random/Randomize.h>
62 #include <CLHEP/Units/GlobalSystemOfUnits.h>
63 #include <CLHEP/Units/GlobalPhysicalConstants.h>
64 
65 #include <cstdint>
66 
67 //#define EDM_ML_DEBUG
68 
69 namespace CLHEP {
70  class HepRandomEngine;
71 }
72 
74  public Observer<const BeginOfRun*>,
75  public Observer<const BeginOfEvent*>,
76  public Observer<const EndOfEvent*>,
77  public Observer<const G4Step*> {
78 public:
80  HcalTB04Analysis(const HcalTB04Analysis&) = delete; // stop default
81  const HcalTB04Analysis& operator=(const HcalTB04Analysis&) = delete;
82  ~HcalTB04Analysis() override;
83 
84  void produce(edm::Event&, const edm::EventSetup&) override;
85 
86 private:
87  void init();
88 
89  // observer methods
90  void update(const BeginOfRun* run) override;
91  void update(const BeginOfEvent* evt) override;
92  void update(const G4Step* step) override;
93  void update(const EndOfEvent* evt) override;
94 
95  //User methods
96  void fillBuffer(const EndOfEvent* evt);
97  void qieAnalysis(CLHEP::HepRandomEngine*);
98  void xtalAnalysis(CLHEP::HepRandomEngine*);
99  void finalAnalysis();
100  void fillEvent(PHcalTB04Info&);
101 
102  void clear();
103  int unitID(uint32_t id);
104  double scale(int det, int layer);
105  double timeOfFlight(int det, int layer, double eta);
106 
107 private:
110 
111  // to read from parameter set
112  bool hcalOnly;
113  int mode, type;
115  int iceta, icphi;
117  std::vector<std::string> names;
118  G4RotationMatrix* beamline_RM;
119 
120  // Constants for the run
121  int count;
123  std::vector<int> idHcal, idXtal;
124  std::vector<uint32_t> idTower, idEcal;
125 
126  // Constants for the event
129  std::vector<CaloHit> ecalHitCache;
130  std::vector<CaloHit> hcalHitCache, hcalHitLayer;
131  std::vector<double> esimh, eqie, esime, enois;
132  std::vector<double> eseta, eqeta, esphi, eqphi, eslay, eqlay;
134 
135  bool pvFound;
136  int evNum, pvType;
137  G4ThreeVector pvPosition, pvMomentum, pvUVW;
138  std::vector<int> secTrackID, secPartID;
139  std::vector<G4ThreeVector> secMomentum;
140  std::vector<double> secEkin;
141  std::vector<int> shortLivedSecondaries;
142 };
143 
144 //
145 // constructors and destructor
146 //
147 
148 HcalTB04Analysis::HcalTB04Analysis(const edm::ParameterSet& p) : myQie(nullptr), histo(nullptr) {
149  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTB04Analysis");
150  hcalOnly = m_Anal.getParameter<bool>("HcalOnly");
151  mode = m_Anal.getParameter<int>("Mode");
152  type = m_Anal.getParameter<int>("Type");
153  ecalNoise = m_Anal.getParameter<double>("EcalNoise");
154  scaleHB0 = m_Anal.getParameter<double>("ScaleHB0");
155  scaleHB16 = m_Anal.getParameter<double>("ScaleHB16");
156  scaleHO = m_Anal.getParameter<double>("ScaleHO");
157  scaleHE0 = m_Anal.getParameter<double>("ScaleHE0");
158  names = m_Anal.getParameter<std::vector<std::string> >("Names");
159  beamOffset = -m_Anal.getParameter<double>("BeamPosition") * cm;
160  double fMinEta = m_Anal.getParameter<double>("MinEta");
161  double fMaxEta = m_Anal.getParameter<double>("MaxEta");
162  double fMinPhi = m_Anal.getParameter<double>("MinPhi");
163  double fMaxPhi = m_Anal.getParameter<double>("MaxPhi");
164  double beamEta = (fMaxEta + fMinEta) / 2.;
165  double beamPhi = (fMaxPhi + fMinPhi) / 2.;
166  double beamThet = 2 * atan(exp(-beamEta));
167  if (beamPhi < 0)
168  beamPhi += twopi;
169  iceta = (int)(beamEta / 0.087) + 1;
170  icphi = (int)(fabs(beamPhi) / 0.087) + 5;
171  if (icphi > 72)
172  icphi -= 73;
173 
174  produces<PHcalTB04Info>();
175 
176  beamline_RM = new G4RotationMatrix;
177  beamline_RM->rotateZ(-beamPhi);
178  beamline_RM->rotateY(-beamThet);
179 
180  edm::LogVerbatim("HcalTBSim")
181  << "HcalTB04:: Initialised as observer of BeginOf Job/BeginOfRun/BeginOfEvent/G4Step/EndOfEvent with Parameter "
182  "values:\n \thcalOnly = "
183  << hcalOnly << "\tecalNoise = " << ecalNoise << "\n\tMode = " << mode << " (0: HB2 Standard; 1:HB2 Segmented)"
184  << "\tType = " << type << " (0: HB; 1 HE; 2 HB+HE)\n\tbeamOffset = " << beamOffset << "\ticeta = " << iceta
185  << "\ticphi = " << icphi << "\n\tbeamline_RM = " << *beamline_RM;
186 
187  init();
188 
189  myQie = new HcalQie(p);
190  histo = new HcalTB04Histo(m_Anal);
191 }
192 
194 #ifdef EDM_ML_DEBUG
195  edm::LogVerbatim("HcalTBSim") << "\n --------> Total number of selected entries : " << count << "\nPointers:: QIE "
196  << myQie << " Histo " << histo;
197 #endif
198  if (myQie) {
199  delete myQie;
200  myQie = nullptr;
201  }
202  if (histo) {
203  delete histo;
204  histo = nullptr;
205  }
206 }
207 
208 //
209 // member functions
210 //
211 
213  std::unique_ptr<PHcalTB04Info> product(new PHcalTB04Info);
214  fillEvent(*product);
215  e.put(std::move(product));
216 }
217 
220  nTower = idTower.size();
221 #ifdef EDM_ML_DEBUG
222  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Save information from " << nTower << " HCal towers";
223 #endif
224  idHcal.reserve(nTower);
225  for (int i = 0; i < nTower; i++) {
226  int id = unitID(idTower[i]);
227  idHcal.push_back(id);
228 #ifdef EDM_ML_DEBUG
229  edm::LogVerbatim("HcalTBSim") << "\tTower[" << i << "] Original " << std::hex << idTower[i] << " Stored "
230  << idHcal[i] << std::dec;
231 #endif
232  }
233 
234  if (!hcalOnly) {
235  int det = 10;
236  uint32_t id1;
237  nCrystal = 0;
238  for (int lay = 1; lay < 8; lay++) {
239  for (int icr = 1; icr < 8; icr++) {
240  id1 = HcalTestNumbering::packHcalIndex(det, 0, 1, icr, lay, 1);
241  int id = unitID(id1);
242  idEcal.push_back(id1);
243  idXtal.push_back(id);
244  nCrystal++;
245  }
246  }
247  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Save information from " << nCrystal << " ECal Crystals";
248 #ifdef EDM_ML_DEBUG
249  for (int i = 0; i < nCrystal; i++) {
250  edm::LogVerbatim("HcalTBSim") << "\tCrystal[" << i << "] Original " << std::hex << idEcal[i] << " Stored "
251  << idXtal[i] << std::dec;
252  }
253 #endif
254  }
255  // Profile vectors
256  eseta.reserve(5);
257  eqeta.reserve(5);
258  esphi.reserve(3);
259  eqphi.reserve(3);
260  eslay.reserve(20);
261  eqlay.reserve(20);
262  for (int i = 0; i < 5; i++) {
263  eseta.push_back(0.);
264  eqeta.push_back(0.);
265  }
266  for (int i = 0; i < 3; i++) {
267  esphi.push_back(0.);
268  eqphi.push_back(0.);
269  }
270  for (int i = 0; i < 20; i++) {
271  eslay.push_back(0.);
272  eqlay.push_back(0.);
273  }
274 
275  // counter
276  count = 0;
277  evNum = 0;
278  clear();
279 }
280 
282  int irun = (*run)()->GetRunID();
283  edm::LogVerbatim("HcalTBSim") << " =====> Begin of Run = " << irun;
284 
285  G4SDManager* sd = G4SDManager::GetSDMpointerIfExist();
286  if (sd != nullptr) {
287  std::string sdname = names[0];
288  G4VSensitiveDetector* aSD = sd->FindSensitiveDetector(sdname);
289  if (aSD == nullptr) {
290  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::beginOfRun: No SD"
291  << " with name " << sdname << " in this "
292  << "Setup";
293  } else {
294  HCalSD* theCaloSD = dynamic_cast<HCalSD*>(aSD);
295  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::beginOfRun: Finds SD with name " << theCaloSD->GetName()
296  << " in this Setup";
298  theCaloSD->setNumberingScheme(org);
299  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::beginOfRun: set a new numbering scheme";
300  }
301  if (!hcalOnly) {
302  sdname = names[1];
303  aSD = sd->FindSensitiveDetector(sdname);
304  if (aSD == nullptr) {
305  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::beginOfRun: No SD"
306  << " with name " << sdname << " in this "
307  << "Setup";
308  } else {
309  ECalSD* theCaloSD = dynamic_cast<ECalSD*>(aSD);
310  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::beginOfRun: Finds SD with name " << theCaloSD->GetName()
311  << " in this Setup";
313  theCaloSD->setNumberingScheme(org);
314  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::beginOfRun: set a new numbering scheme";
315  }
316  }
317  } else {
318  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::beginOfRun: Could "
319  << "not get SD Manager!";
320  }
321 }
322 
324  clear();
325 #ifdef EDM_ML_DEBUG
326  evNum = (*evt)()->GetEventID();
327  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis: =====> Begin of event = " << evNum;
328 #endif
329 }
330 
331 void HcalTB04Analysis::update(const G4Step* aStep) {
332  if (aStep != nullptr) {
333  //Get Step properties
334  G4ThreeVector thePreStepPoint = aStep->GetPreStepPoint()->GetPosition();
335  G4ThreeVector thePostStepPoint;
336 
337  // Get Tracks properties
338  G4Track* aTrack = aStep->GetTrack();
339  int trackID = aTrack->GetTrackID();
340  int parentID = aTrack->GetParentID();
341  const G4ThreeVector& position = aTrack->GetPosition();
342  G4ThreeVector momentum = aTrack->GetMomentum();
343  G4String partType = aTrack->GetDefinition()->GetParticleType();
344  G4String partSubType = aTrack->GetDefinition()->GetParticleSubType();
345  int partPDGEncoding = aTrack->GetDefinition()->GetPDGEncoding();
346 #ifdef EDM_ML_DEBUG
347  bool isPDGStable = aTrack->GetDefinition()->GetPDGStable();
348 #endif
349  double pDGlifetime = aTrack->GetDefinition()->GetPDGLifeTime();
350  double gammaFactor = aStep->GetPreStepPoint()->GetGamma();
351 
352  if (!pvFound) { //search for v1
353  double stepDeltaEnergy = aStep->GetDeltaEnergy();
354  double kinEnergy = aTrack->GetKineticEnergy();
355 
356  // look for DeltaE > 10% kinEnergy of particle, or particle death - Ek=0
357  if (trackID == 1 && parentID == 0 && ((kinEnergy == 0.) || (fabs(stepDeltaEnergy / kinEnergy) > 0.1))) {
358  pvType = -1;
359  if (kinEnergy == 0.) {
360  pvType = 0;
361  } else {
362  if (fabs(stepDeltaEnergy / kinEnergy) > 0.1)
363  pvType = 1;
364  }
365  pvFound = true;
367  pvMomentum = momentum;
368  // Rotated coord.system:
369  pvUVW = (*beamline_RM) * (pvPosition);
370 
371  //Volume name requires some checks:
372  G4String thePostPVname = "NoName";
373  G4StepPoint* thePostPoint = aStep->GetPostStepPoint();
374  if (thePostPoint) {
375  thePostStepPoint = thePostPoint->GetPosition();
376  G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
377  if (thePostPV)
378  thePostPVname = thePostPV->GetName();
379  }
380 #ifdef EDM_ML_DEBUG
381  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: V1 found at: " << thePostStepPoint
382  << " G4VPhysicalVolume: " << thePostPVname;
383  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::fill_v1Pos: Primary Track momentum: " << pvMomentum
384  << " psoition " << pvPosition << " u/v/w " << pvUVW;
385 #endif
386  }
387  } else {
388  // watch for secondaries originating @v1, including the surviving primary
389  if ((trackID != 1 && parentID == 1 && (aTrack->GetCurrentStepNumber() == 1) && (thePreStepPoint == pvPosition)) ||
390  (trackID == 1 && thePreStepPoint == pvPosition)) {
391 #ifdef EDM_ML_DEBUG
392  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::A secondary... PDG:" << partPDGEncoding
393  << " TrackID:" << trackID << " ParentID:" << parentID
394  << " stable: " << isPDGStable << " Tau: " << pDGlifetime
395  << " cTauGamma=" << c_light * pDGlifetime * gammaFactor * 1000.
396  << "um GammaFactor: " << gammaFactor;
397 #endif
398  secTrackID.push_back(trackID);
399  secPartID.push_back(partPDGEncoding);
400  secMomentum.push_back(momentum);
401  secEkin.push_back(aTrack->GetKineticEnergy());
402 
403  // Check for short-lived secondaries: cTauGamma<100um
404  double ctaugamma_um = CLHEP::c_light * pDGlifetime * gammaFactor * 1000.;
405  if ((ctaugamma_um > 0.) && (ctaugamma_um < 100.)) { //short-lived secondary
406  shortLivedSecondaries.push_back(trackID);
407  } else { //normal secondary - enter into the V1-calorimetric tree
408  // histos->fill_v1cSec (aTrack);
409  }
410  }
411  // Also watch for tertiary particles coming from
412  // short-lived secondaries from V1
413  if (aTrack->GetCurrentStepNumber() == 1) {
414  if (!shortLivedSecondaries.empty()) {
415  int pid = parentID;
416  std::vector<int>::iterator pos1 = shortLivedSecondaries.begin();
417  std::vector<int>::iterator pos2 = shortLivedSecondaries.end();
418  std::vector<int>::iterator pos;
419  for (pos = pos1; pos != pos2; pos++) {
420  if (*pos == pid) { //ParentID is on the list of short-lived
421  // secondary
422 #ifdef EDM_ML_DEBUG
423  edm::LogVerbatim("HcalTBSim")
424  << "HcalTB04Analysis:: A tertiary... PDG:" << partPDGEncoding << " TrackID:" << trackID
425  << " ParentID:" << parentID << " stable: " << isPDGStable << " Tau: " << pDGlifetime
426  << " cTauGamma=" << c_light * pDGlifetime * gammaFactor * 1000. << "um GammaFactor: " << gammaFactor;
427 #endif
428  }
429  }
430  }
431  }
432  }
433  }
434 }
435 
437  count++;
438 
439  //fill the buffer
440 #ifdef EDM_ML_DEBUG
441  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::Fill event " << (*evt)()->GetEventID();
442 #endif
443  fillBuffer(evt);
444 
445  //QIE analysis
446 #ifdef EDM_ML_DEBUG
447  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::Do QIE analysis with " << hcalHitCache.size() << " hits";
448 #endif
449  CLHEP::HepRandomEngine* engine = G4Random::getTheEngine();
450  qieAnalysis(engine);
451 
452  //Energy in Crystal Matrix
453  if (!hcalOnly) {
454 #ifdef EDM_ML_DEBUG
455  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::Do Xtal analysis with " << ecalHitCache.size() << " hits";
456 #endif
457  xtalAnalysis(engine);
458  }
459 
460  //Final Analysis
461 #ifdef EDM_ML_DEBUG
462  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::Final analysis";
463 #endif
464  finalAnalysis();
465 
466  int iEvt = (*evt)()->GetEventID();
467  if (iEvt < 10)
468  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
469  else if ((iEvt < 100) && (iEvt % 10 == 0))
470  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
471  else if ((iEvt < 1000) && (iEvt % 100 == 0))
472  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
473  else if ((iEvt < 10000) && (iEvt % 1000 == 0))
474  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
475 }
476 
478  std::vector<CaloHit> hhits, hhitl;
479  int idHC, j;
480  CaloG4HitCollection* theHC;
481  std::map<int, float, std::less<int> > primaries;
482  double etot1 = 0, etot2 = 0;
483 
484  // Look for the Hit Collection of HCal
485  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
486  std::string sdName = names[0];
487  idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
488  theHC = (CaloG4HitCollection*)allHC->GetHC(idHC);
489 #ifdef EDM_ML_DEBUG
490  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Hit Collection for " << sdName << " of ID " << idHC
491  << " is obtained at " << theHC << " with " << theHC->entries() << " entries";
492 #endif
493  int thehc_entries = theHC->entries();
494  if (idHC >= 0 && theHC != nullptr) {
495  hhits.reserve(theHC->entries());
496  hhitl.reserve(theHC->entries());
497  for (j = 0; j < thehc_entries; j++) {
498  CaloG4Hit* aHit = (*theHC)[j];
499  double e = aHit->getEnergyDeposit() / CLHEP::GeV;
500  double time = aHit->getTimeSlice();
501  math::XYZPoint pos = aHit->getEntry();
502  unsigned int id = aHit->getUnitID();
503  double theta = pos.theta();
504  double eta = -log(tan(theta * 0.5));
505  double phi = pos.phi();
506  int det, z, group, ieta, iphi, layer;
507  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, layer);
508  double jitter = time - timeOfFlight(det, layer, eta);
509  if (jitter < 0)
510  jitter = 0;
511  if (e < 0 || e > 1.)
512  e = 0;
513  double escl = e * scale(det, layer);
514  unsigned int idx = HcalTBNumberingScheme::getUnitID(id, mode);
515  CaloHit hit(det, layer, escl, eta, phi, jitter, idx);
516  hhits.push_back(hit);
517  CaloHit hitl(det, layer, escl, eta, phi, jitter, id);
518  hhitl.push_back(hitl);
519  primaries[aHit->getTrackID()] += e;
520  etot1 += escl;
521 #ifdef EDM_ML_DEBUG
522  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Hcal Hit i/p " << j << " ID 0x" << std::hex << id << " 0x"
523  << idx << std::dec << " time " << std::setw(6) << time << " " << std::setw(6)
524  << jitter << " theta " << std::setw(8) << theta << " eta " << std::setw(8) << eta
525  << " phi " << std::setw(8) << phi << " e " << std::setw(8) << e << " "
526  << std::setw(8) << escl;
527 #endif
528  }
529  }
530 
531  // Add hits in the same channel within same time slice
532  std::vector<CaloHit>::iterator itr;
533  int nHit = hhits.size();
534  std::vector<CaloHit*> hits(nHit);
535  for (j = 0, itr = hhits.begin(); itr != hhits.end(); j++, itr++) {
536  hits[j] = &hhits[j];
537  }
538  sort(hits.begin(), hits.end(), CaloHitIdMore());
539  std::vector<CaloHit*>::iterator k1, k2;
540  int nhit = 0;
541  for (k1 = hits.begin(); k1 != hits.end(); k1++) {
542  int det = (**k1).det();
543  int layer = (**k1).layer();
544  double ehit = (**k1).e();
545  double eta = (**k1).eta();
546  double phi = (**k1).phi();
547  double jitter = (**k1).t();
548  uint32_t unitID = (**k1).id();
549  int jump = 0;
550  for (k2 = k1 + 1; k2 != hits.end() && fabs(jitter - (**k2).t()) < 1 && unitID == (**k2).id(); k2++) {
551  ehit += (**k2).e();
552  jump++;
553  }
554  nhit++;
555  CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
556  hcalHitCache.push_back(hit);
557  etot2 += ehit;
558  k1 += jump;
559 #ifdef EDM_ML_DEBUG
560  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Hcal Hit store " << nhit << " ID 0x" << std::hex << unitID
561  << std::dec << " time " << std::setw(6) << jitter << " eta " << std::setw(8) << eta
562  << " phi " << std::setw(8) << phi << " e " << std::setw(8) << ehit;
563 #endif
564  }
565 #ifdef EDM_ML_DEBUG
566  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Stores " << nhit << " HCal hits from " << nHit
567  << " input hits E(Hcal) " << etot1 << " " << etot2;
568 #endif
569  //Repeat for Hit in each layer (hhits and hhitl sizes are the same)
570  for (j = 0, itr = hhitl.begin(); itr != hhitl.end(); j++, itr++) {
571  hits[j] = &hhitl[j];
572  }
573  sort(hits.begin(), hits.end(), CaloHitIdMore());
574  int nhitl = 0;
575  double etotl = 0;
576  for (k1 = hits.begin(); k1 != hits.end(); k1++) {
577  int det = (**k1).det();
578  int layer = (**k1).layer();
579  double ehit = (**k1).e();
580  double eta = (**k1).eta();
581  double phi = (**k1).phi();
582  double jitter = (**k1).t();
583  uint32_t unitID = (**k1).id();
584  int jump = 0;
585  for (k2 = k1 + 1; k2 != hits.end() && fabs(jitter - (**k2).t()) < 1 && unitID == (**k2).id(); k2++) {
586  ehit += (**k2).e();
587  jump++;
588  }
589  nhitl++;
590  CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
591  hcalHitLayer.push_back(hit);
592  etotl += ehit;
593  k1 += jump;
594 #ifdef EDM_ML_DEBUG
595  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Hcal Hit store " << nhitl << " ID 0x" << std::hex << unitID
596  << std::dec << " time " << std::setw(6) << jitter << " eta " << std::setw(8) << eta
597  << " phi " << std::setw(8) << phi << " e " << std::setw(8) << ehit;
598 #endif
599  }
600 #ifdef EDM_ML_DEBUG
601  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Stores " << nhitl << " HCal hits from " << nHit
602  << " input hits E(Hcal) " << etot1 << " " << etotl;
603 #endif
604  // Look for the Hit Collection of ECal
605  std::vector<CaloHit> ehits;
606  sdName = names[1];
607  idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
608  theHC = (CaloG4HitCollection*)allHC->GetHC(idHC);
609  etot1 = etot2 = 0;
610 #ifdef EDM_ML_DEBUG
611  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Hit Collection for " << sdName << " of ID " << idHC
612  << " is obtained at " << theHC << " with " << theHC->entries() << " entries";
613 #endif
614  if (idHC >= 0 && theHC != nullptr) {
615  thehc_entries = theHC->entries();
616  ehits.reserve(theHC->entries());
617  for (j = 0; j < thehc_entries; j++) {
618  CaloG4Hit* aHit = (*theHC)[j];
619  double e = aHit->getEnergyDeposit() / CLHEP::GeV;
620  double time = aHit->getTimeSlice();
621  if (e < 0 || e > 100000.)
622  e = 0;
623  if (e > 0) {
624  math::XYZPoint pos = aHit->getEntry();
625  unsigned int id = aHit->getUnitID();
626  double theta = pos.theta();
627  double eta = -log(tan(theta * 0.5));
628  double phi = pos.phi();
629  int det, z, group, ieta, iphi, layer;
630  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, layer);
631  CaloHit hit(det, 0, e, eta, phi, time, id);
632  ehits.push_back(hit);
633  primaries[aHit->getTrackID()] += e;
634  etot1 += e;
635 #ifdef EDM_ML_DEBUG
636  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Ecal Hit i/p " << j << " ID 0x" << std::hex << id
637  << std::dec << " time " << std::setw(6) << time << " theta " << std::setw(8)
638  << theta << " eta " << std::setw(8) << eta << " phi " << std::setw(8) << phi
639  << " e " << std::setw(8) << e;
640 #endif
641  }
642  }
643  }
644 
645  // Add hits in the same channel within same time slice
646  nHit = ehits.size();
647  std::vector<CaloHit*> hite(nHit);
648  for (j = 0, itr = ehits.begin(); itr != ehits.end(); j++, itr++) {
649  hite[j] = &ehits[j];
650  }
651  sort(hite.begin(), hite.end(), CaloHitIdMore());
652  nhit = 0;
653  for (k1 = hite.begin(); k1 != hite.end(); k1++) {
654  int det = (**k1).det();
655  int layer = (**k1).layer();
656  double ehit = (**k1).e();
657  double eta = (**k1).eta();
658  double phi = (**k1).phi();
659  double jitter = (**k1).t();
660  uint32_t unitID = (**k1).id();
661  int jump = 0;
662  for (k2 = k1 + 1; k2 != hite.end() && fabs(jitter - (**k2).t()) < 1 && unitID == (**k2).id(); k2++) {
663  ehit += (**k2).e();
664  jump++;
665  }
666  nhit++;
667  CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
668  ecalHitCache.push_back(hit);
669  etot2 += ehit;
670  k1 += jump;
671 #ifdef EDM_ML_DEBUG
672  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Ecal Hit store " << nhit << " ID 0x" << std::hex << unitID
673  << std::dec << " time " << std::setw(6) << jitter << " eta " << std::setw(8) << eta
674  << " phi " << std::setw(8) << phi << " e " << std::setw(8) << ehit;
675 #endif
676  }
677 #ifdef EDM_ML_DEBUG
678  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Stores " << nhit << " ECal hits from " << nHit
679  << " input hits E(Ecal) " << etot1 << " " << etot2;
680 #endif
681  // Find Primary info:
682  nPrimary = (int)(primaries.size());
683  int trackID = 0;
684  G4PrimaryParticle* thePrim = nullptr;
685  int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
686 #ifdef EDM_ML_DEBUG
687  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Event has " << nvertex << " verteices";
688 #endif
689  if (nvertex <= 0)
690  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::EndOfEvent ERROR: no vertex found for event " << evNum;
691  for (int i = 0; i < nvertex; i++) {
692  G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
693  if (avertex == nullptr) {
694  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::EndOfEvent ERR: pointer to vertex = 0 for event " << evNum;
695  } else {
696  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::Vertex number :" << i << " " << avertex->GetPosition();
697  int npart = avertex->GetNumberOfParticle();
698  if (npart == 0)
699  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::End Of Event ERR: no primary!";
700  if (thePrim == nullptr)
701  thePrim = avertex->GetPrimary(trackID);
702  }
703  }
704 
705  if (thePrim != nullptr) {
706  double px = thePrim->GetPx();
707  double py = thePrim->GetPy();
708  double pz = thePrim->GetPz();
709  double p = std::sqrt(pow(px, 2.) + pow(py, 2.) + pow(pz, 2.));
710  pInit = p / CLHEP::GeV;
711  if (p == 0)
712  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis:: EndOfEvent ERR: primary has p=0 ";
713  else {
714  double costheta = pz / p;
715  double theta = acos(std::min(std::max(costheta, -1.), 1.));
716  etaInit = -log(tan(theta / 2));
717  if (px != 0 || py != 0)
718  phiInit = atan2(py, px);
719  }
720  particleType = thePrim->GetPDGcode();
721  } else
722  edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::EndOfEvent ERR: could not find primary";
723 }
724 
725 void HcalTB04Analysis::qieAnalysis(CLHEP::HepRandomEngine* engine) {
726  int hittot = hcalHitCache.size();
727  if (hittot <= 0)
728  hittot = 1;
729  std::vector<CaloHit> hits(hittot);
730  std::vector<int> todo(nTower, 0);
731 
732 #ifdef EDM_ML_DEBUG
733  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::qieAnalysis: Size " << hits.size() << " " << todo.size() << " "
734  << idTower.size() << " " << esimh.size() << " " << eqie.size();
735 #endif
736  // Loop over all HCal hits
737  for (unsigned int k1 = 0; k1 < hcalHitCache.size(); k1++) {
738  CaloHit hit = hcalHitCache[k1];
739  uint32_t id = hit.id();
740  int nhit = 0;
741  double esim = hit.e();
742  hits[nhit] = hit;
743  for (unsigned int k2 = k1 + 1; k2 < hcalHitCache.size(); k2++) {
744  hit = hcalHitCache[k2];
745  if (hit.id() == id) {
746  nhit++;
747  hits[nhit] = hit;
748  esim += hit.e();
749  }
750  }
751  k1 += nhit;
752  nhit++;
753  std::vector<int> cd = myQie->getCode(nhit, hits, engine);
754  double eq = myQie->getEnergy(cd);
755 #ifdef EDM_ML_DEBUG
756  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex << id << std::dec << " registers " << esim
757  << " energy from " << nhit << " hits starting with hit # " << k1
758  << " energy with noise " << eq;
759 #endif
760  for (int k2 = 0; k2 < nTower; k2++) {
761  if (id == idTower[k2]) {
762  todo[k2] = 1;
763  esimh[k2] = esim;
764  eqie[k2] = eq;
765  }
766  }
767  }
768 
769  // Towers with no hit
770  for (int k2 = 0; k2 < nTower; k2++) {
771  if (todo[k2] == 0) {
772  std::vector<int> cd = myQie->getCode(0, hits, engine);
773  double eq = myQie->getEnergy(cd);
774  esimh[k2] = 0;
775  eqie[k2] = eq;
776 #ifdef EDM_ML_DEBUG
777  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex << idTower[k2] << std::dec
778  << " registers " << esimh[k2] << " energy from hits and energy after QIE analysis "
779  << eqie[k2];
780 #endif
781  }
782  }
783 }
784 
785 void HcalTB04Analysis::xtalAnalysis(CLHEP::HepRandomEngine* engine) {
786  CLHEP::RandGaussQ randGauss(*engine);
787 
788  // Crystal Data
789  std::vector<int> iok(nCrystal, 0);
790 #ifdef EDM_ML_DEBUG
791  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::xtalAnalysis: Size " << iok.size() << " " << idEcal.size() << " "
792  << esime.size() << " " << enois.size();
793 #endif
794  for (unsigned int k1 = 0; k1 < ecalHitCache.size(); k1++) {
795  uint32_t id = ecalHitCache[k1].id();
796  int nhit = 0;
797  double esim = ecalHitCache[k1].e();
798  for (unsigned int k2 = k1 + 1; k2 < ecalHitCache.size(); k2++) {
799  if (ecalHitCache[k2].id() == id) {
800  nhit++;
801  esim += ecalHitCache[k2].e();
802  }
803  }
804  k1 += nhit;
805  nhit++;
806  double eq = esim + randGauss.fire(0., ecalNoise);
807 #ifdef EDM_ML_DEBUG
808  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex << id << std::dec << " registers " << esim
809  << " energy from " << nhit << " hits starting with hit # " << k1
810  << " energy with noise " << eq;
811 #endif
812  for (int k2 = 0; k2 < nCrystal; k2++) {
813  if (id == idEcal[k2]) {
814  iok[k2] = 1;
815  esime[k2] = esim;
816  enois[k2] = eq;
817  }
818  }
819  }
820 
821  // Crystals with no hit
822  for (int k2 = 0; k2 < nCrystal; k2++) {
823  if (iok[k2] == 0) {
824  esime[k2] = 0;
825  enois[k2] = randGauss.fire(0., ecalNoise);
826 #ifdef EDM_ML_DEBUG
827  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex << idEcal[k2] << std::dec
828  << " registers " << esime[k2] << " energy from hits and energy from noise "
829  << enois[k2];
830 #endif
831  }
832  }
833 }
834 
836  //Beam Information
838 
839  // Total Energy
840  eecals = ehcals = eecalq = ehcalq = 0.;
841  for (int i = 0; i < nTower; i++) {
842  ehcals += esimh[i];
843  ehcalq += eqie[i];
844  }
845  for (int i = 0; i < nCrystal; i++) {
846  eecals += esime[i];
847  eecalq += enois[i];
848  }
849  etots = eecals + ehcals;
850  etotq = eecalq + ehcalq;
851 #ifdef EDM_ML_DEBUG
852  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Energy deposit at Sim Level (Total) " << etots << " (ECal) "
853  << eecals << " (HCal) " << ehcals
854  << "\nHcalTB04Analysis:: Energy deposit at Qie Level (Total) " << etotq << " (ECal) "
855  << eecalq << " (HCal) " << ehcalq;
856 #endif
857  histo->fillEdep(etots, eecals, ehcals, etotq, eecalq, ehcalq);
858 
859  // Lateral Profile
860  for (int i = 0; i < 5; i++) {
861  eseta[i] = 0.;
862  eqeta[i] = 0.;
863  }
864  for (int i = 0; i < 3; i++) {
865  esphi[i] = 0.;
866  eqphi[i] = 0.;
867  }
868  double e1 = 0, e2 = 0;
869  unsigned int id;
870  for (int i = 0; i < nTower; i++) {
871  int det, z, group, ieta, iphi, layer;
872  id = idTower[i];
873  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, layer);
874  iphi -= (icphi - 1);
875  if (icphi > 4) {
876  if (ieta == 0)
877  ieta = 2;
878  else
879  ieta = -1;
880  } else {
881  ieta = ieta - iceta + 2;
882  }
883  if (iphi >= 0 && iphi < 3 && ieta >= 0 && ieta < 5) {
884  eseta[ieta] += esimh[i];
885  esphi[iphi] += esimh[i];
886  e1 += esimh[i];
887  eqeta[ieta] += eqie[i];
888  eqphi[iphi] += eqie[i];
889  e2 += eqie[i];
890  }
891  }
892  for (int i = 0; i < 3; i++) {
893  if (e1 > 0)
894  esphi[i] /= e1;
895  if (e2 > 0)
896  eqphi[i] /= e2;
897  }
898  for (int i = 0; i < 5; i++) {
899  if (e1 > 0)
900  eseta[i] /= e1;
901  if (e2 > 0)
902  eqeta[i] /= e2;
903  }
904 #ifdef EDM_ML_DEBUG
905  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Energy fraction along Eta and Phi (Sim/Qie)";
906  for (int i = 0; i < 5; i++)
907  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: [" << i << "] Eta Sim = " << eseta[i] << " Qie = " << eqeta[i]
908  << " Phi Sim = " << esphi[i] << " Qie = " << eqphi[i];
909 #endif
910  histo->fillTrnsProf(eseta, eqeta, esphi, eqphi);
911 
912  // Longitudianl profile
913  for (int i = 0; i < 20; i++) {
914  eslay[i] = 0.;
915  eqlay[i] = 0.;
916  }
917  e1 = 0;
918  e2 = 0;
919  for (int i = 0; i < nTower; i++) {
920  int det, z, group, ieta, iphi, layer;
921  id = idTower[i];
922  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, layer);
923  iphi -= (icphi - 1);
924  layer -= 1;
925  if (iphi >= 0 && iphi < 3 && layer >= 0 && layer < 20) {
926  eslay[layer] += esimh[i];
927  e1 += esimh[i];
928  eqlay[layer] += eqie[i];
929  e2 += eqie[i];
930  }
931  }
932  for (int i = 0; i < 20; i++) {
933  if (e1 > 0)
934  eslay[i] /= e1;
935  if (e2 > 0)
936  eqlay[i] /= e2;
937  }
938 #ifdef EDM_ML_DEBUG
939  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Energy fraction along Layer";
940  for (int i = 0; i < 20; i++)
941  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: [" << i << "] Sim = " << eslay[i] << " Qie = " << eqlay[i];
942 #endif
943  histo->fillLongProf(eslay, eqlay);
944 }
945 
947  //Setup the ID's
948  product.setIDs(idHcal, idXtal);
949 
950  //Beam Information
952 
953  //Energy deposits in the crystals and towers
954  product.setEdepHcal(esimh, eqie);
955  product.setEdepHcal(esime, enois);
956 
957  // Total Energy
958  product.setEdep(etots, eecals, ehcals, etotq, eecalq, ehcalq);
959 
960  // Lateral Profile
961  product.setTrnsProf(eseta, eqeta, esphi, eqphi);
962 
963  // Longitudianl profile
964  product.setLongProf(eslay, eqlay);
965 
966  //Save Hits
967  int i, nhit = 0;
968  std::vector<CaloHit>::iterator itr;
969  for (i = 0, itr = ecalHitCache.begin(); itr != ecalHitCache.end(); i++, itr++) {
970  uint32_t id = itr->id();
971  int det, z, group, ieta, iphi, lay;
972  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, lay);
973  product.saveHit(det, lay, ieta, iphi, itr->e(), itr->t());
974  nhit++;
975 #ifdef EDM_ML_DEBUG
976  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Save Hit " << std::setw(3) << i + 1 << " ID 0x" << std::hex
977  << group << std::dec << " " << std::setw(2) << det << " " << std::setw(2) << lay
978  << " " << std::setw(1) << z << " " << std::setw(3) << ieta << " " << std::setw(3)
979  << iphi << " T " << std::setw(6) << itr->t() << " E " << std::setw(6) << itr->e();
980 #endif
981  }
982  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Saves " << nhit << " hits from Crystals";
983  int hit = nhit;
984  nhit = 0;
985 
986  for (i = hit, itr = hcalHitCache.begin(); itr != hcalHitCache.end(); i++, itr++) {
987  uint32_t id = itr->id();
988  int det, z, group, ieta, iphi, lay;
989  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, lay);
990  product.saveHit(det, lay, ieta, iphi, itr->e(), itr->t());
991  nhit++;
992 #ifdef EDM_ML_DEBUG
993  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Save Hit " << std::setw(3) << i + 1 << " ID 0x" << std::hex
994  << group << std::dec << " " << std::setw(2) << det << " " << std::setw(2) << lay
995  << " " << std::setw(1) << z << " " << std::setw(3) << ieta << " " << std::setw(3)
996  << iphi << " T " << std::setw(6) << itr->t() << " E " << std::setw(6) << itr->e();
997 #endif
998  }
999  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis:: Saves " << nhit << " hits from HCal";
1000 
1001  //Vertex associated quantities
1002  product.setVtxPrim(evNum,
1003  pvType,
1004  pvPosition.x(),
1005  pvPosition.y(),
1006  pvPosition.z(),
1007  pvUVW.x(),
1008  pvUVW.y(),
1009  pvUVW.z(),
1010  pvMomentum.x(),
1011  pvMomentum.y(),
1012  pvMomentum.z());
1013  for (unsigned int i = 0; i < secTrackID.size(); i++) {
1014  product.setVtxSec(
1015  secTrackID[i], secPartID[i], secMomentum[i].x(), secMomentum[i].y(), secMomentum[i].z(), secEkin[i]);
1016  }
1017 }
1018 
1020  pvFound = false;
1021  pvType = -2;
1022  pvPosition = G4ThreeVector();
1023  pvMomentum = G4ThreeVector();
1024  pvUVW = G4ThreeVector();
1025  secTrackID.clear();
1026  secPartID.clear();
1027  secMomentum.clear();
1028  secEkin.clear();
1029  shortLivedSecondaries.clear();
1030 
1031  ecalHitCache.erase(ecalHitCache.begin(), ecalHitCache.end());
1032  hcalHitCache.erase(hcalHitCache.begin(), hcalHitCache.end());
1033  hcalHitLayer.erase(hcalHitLayer.begin(), hcalHitLayer.end());
1034  nPrimary = particleType = 0;
1035  pInit = etaInit = phiInit = 0;
1036 
1037  esimh.clear();
1038  eqie.clear();
1039  esimh.reserve(nTower);
1040  eqie.reserve(nTower);
1041  for (int i = 0; i < nTower; i++) {
1042  esimh.push_back(0.);
1043  eqie.push_back(0.);
1044  }
1045  esime.clear();
1046  enois.clear();
1047  esime.reserve(nCrystal);
1048  enois.reserve(nCrystal);
1049  for (int i = 0; i < nCrystal; i++) {
1050  esime.push_back(0.);
1051  enois.push_back(0.);
1052  }
1053 }
1054 
1055 int HcalTB04Analysis::unitID(uint32_t id) {
1056  int det, z, group, ieta, iphi, lay;
1057  HcalTestNumbering::unpackHcalIndex(id, det, z, group, ieta, iphi, lay);
1058  group = (det & 15) << 20;
1059  group += ((lay - 1) & 31) << 15;
1060  group += (z & 1) << 14;
1061  group += (ieta & 127) << 7;
1062  group += (iphi & 127);
1063  return group;
1064 }
1065 
1066 double HcalTB04Analysis::scale(int det, int layer) {
1067  double tmp = 1.;
1068  if (det == static_cast<int>(HcalBarrel)) {
1069  if (layer == 1)
1070  tmp = scaleHB0;
1071  else if (layer == 17)
1072  tmp = scaleHB16;
1073  else if (layer > 17)
1074  tmp = scaleHO;
1075  } else {
1076  if (layer <= 2)
1077  tmp = scaleHE0;
1078  }
1079  return tmp;
1080 }
1081 
1082 double HcalTB04Analysis::timeOfFlight(int det, int layer, double eta) {
1083  double theta = 2.0 * atan(exp(-eta));
1084  double dist = beamOffset;
1085  if (det == static_cast<int>(HcalBarrel)) {
1086  const double rLay[19] = {1836.0,
1087  1902.0,
1088  1962.0,
1089  2022.0,
1090  2082.0,
1091  2142.0,
1092  2202.0,
1093  2262.0,
1094  2322.0,
1095  2382.0,
1096  2448.0,
1097  2514.0,
1098  2580.0,
1099  2646.0,
1100  2712.0,
1101  2776.0,
1102  2862.5,
1103  3847.0,
1104  4052.0};
1105  if (layer > 0 && layer <= 19)
1106  dist += rLay[layer - 1] * mm / sin(theta);
1107  } else {
1108  const double zLay[19] = {4034.0,
1109  4032.0,
1110  4123.0,
1111  4210.0,
1112  4297.0,
1113  4384.0,
1114  4471.0,
1115  4558.0,
1116  4645.0,
1117  4732.0,
1118  4819.0,
1119  4906.0,
1120  4993.0,
1121  5080.0,
1122  5167.0,
1123  5254.0,
1124  5341.0,
1125  5428.0,
1126  5515.0};
1127  if (layer > 0 && layer <= 19)
1128  dist += zLay[layer - 1] * mm / cos(theta);
1129  }
1130 
1131  double tmp = dist / c_light / ns;
1132 #ifdef EDM_ML_DEBUG
1133  edm::LogVerbatim("HcalTBSim") << "HcalTB04Analysis::timeOfFlight " << tmp << " for det/lay " << det << " " << layer
1134  << " eta/theta " << eta << " " << theta / deg << " dist " << dist;
1135 #endif
1136  return tmp;
1137 }
1138 
Log< level::Info, true > LogVerbatim
#define DEFINE_SIMWATCHER(type)
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::vector< double > secEkin
void setLongProf(const std::vector< double > &es, const std::vector< double > &eq)
uint16_t *__restrict__ id
std::vector< int > idHcal
void fillPrimary(double energy, double eta, double phi)
std::vector< CaloHit > hcalHitLayer
std::vector< double > eqeta
std::vector< double > eqie
void fillTrnsProf(const std::vector< double > &es1, const std::vector< double > &eq1, const std::vector< double > &es2, const std::vector< double > &eq2)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void fillBuffer(const EndOfEvent *evt)
std::vector< int > shortLivedSecondaries
Geom::Theta< T > theta() const
std::vector< int > secTrackID
double npart
Definition: HydjetWrapper.h:46
std::vector< double > eqphi
void setPrimary(int primary, int id, double energy, double eta, double phi)
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:544
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
void setEdep(double simtot, double sime, double simh, double digtot, double dige, double digh)
void setTrnsProf(const std::vector< double > &es1, const std::vector< double > &eq1, const std::vector< double > &es2, const std::vector< double > &eq2)
constexpr std::array< uint8_t, layerIndexSize > layer
void saveHit(int det, int lay, int eta, int phi, double e, double t)
std::vector< int > secPartID
Definition: HCalSD.h:38
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
std::vector< uint32_t > idEcal
std::vector< double > eseta
void setIDs(const std::vector< int > &, const std::vector< int > &)
uint32_t id() const
Definition: CaloHit.h:25
void setVtxSec(int id, int pdg, double px, double py, double pz, double ek)
G4ThreeVector pvUVW
void setEdepHcal(const std::vector< double > &esim, const std::vector< double > &edig)
double getEnergy(const std::vector< int > &)
Definition: HcalQie.cc:354
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)
G4RotationMatrix * beamline_RM
std::vector< double > esime
T min(T a, T b)
Definition: MathUtil.h:58
static std::vector< uint32_t > getUnitIDs(const int type, const int mode)
double scale(int det, int layer)
std::vector< double > esphi
double timeOfFlight(int det, int layer, double eta)
unsigned int id
Definition: ECalSD.h:31
tuple group
Definition: watchdog.py:82
std::vector< CaloHit > hcalHitCache
int getTrackID() const
Definition: CaloG4Hit.h:64
void xtalAnalysis(CLHEP::HepRandomEngine *)
static uint32_t getUnitID(const uint32_t id, const int mode)
HcalTB04Analysis(const edm::ParameterSet &p)
~HcalTB04Analysis() override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double sd
std::vector< int > idXtal
std::vector< CaloHit > ecalHitCache
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double e() const
Definition: CaloHit.h:21
std::vector< double > enois
void setVtxPrim(int evNum, int type, double x, double y, double z, double u, double v, double w, double px, double py, double pz)
std::vector< G4ThreeVector > secMomentum
int unitID(uint32_t id)
HcalTB04Histo * histo
std::vector< int > getCode(int, const std::vector< CaloHit > &, CLHEP::HepRandomEngine *)
Definition: HcalQie.cc:268
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
static int position[264][3]
Definition: ReadPGInfo.cc:289
G4ThreeVector pvMomentum
double getTimeSlice() const
Definition: CaloG4Hit.h:67
std::vector< uint32_t > idTower
void qieAnalysis(CLHEP::HepRandomEngine *)
std::vector< double > eslay
math::XYZPoint getEntry() const
Definition: CaloG4Hit.h:46
std::vector< double > eqlay
step
Definition: StallMonitor.cc:94
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
Log< level::Warning, false > LogWarning
G4ThreeVector pvPosition
void setNumberingScheme(EcalNumberingScheme *)
Definition: ECalSD.cc:327
tmp
align.sh
Definition: createJobs.py:716
std::vector< std::string > names
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
void fillEdep(double etots, double eecals, double ehcals, double etotq, double eecalq, double ehcalq)
void fillEvent(PHcalTB04Info &)
const HcalTB04Analysis & operator=(const HcalTB04Analysis &)=delete
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
std::vector< double > esimh
void fillLongProf(const std::vector< double > &es, const std::vector< double > &eq)
void produce(edm::Event &, const edm::EventSetup &) override