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