CMS 3D CMS Logo

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