CMS 3D CMS Logo

HcalTestAnalysis.cc
Go to the documentation of this file.
6 
7 // to retreive hits
12 
15 
19 
25 
26 #include "G4SDManager.hh"
27 #include "G4Step.hh"
28 #include "G4Track.hh"
29 #include "G4VProcess.hh"
30 #include "G4HCofThisEvent.hh"
31 
32 #include <CLHEP/Random/Randomize.h>
33 #include <CLHEP/Units/GlobalSystemOfUnits.h>
34 #include <CLHEP/Units/GlobalPhysicalConstants.h>
35 
36 #include <cmath>
37 #include <iomanip>
38 #include <iostream>
39 #include <string>
40 #include <vector>
41 
43  public Observer<const BeginOfRun*>,
44  public Observer<const BeginOfEvent*>,
45  public Observer<const EndOfEvent*>,
46  public Observer<const G4Step*> {
47 public:
49  ~HcalTestAnalysis() override;
50 
52  void produce(edm::Event&, const edm::EventSetup&) override;
53  void beginRun(edm::EventSetup const&) override;
54 
55 private:
56  // observer classes
57  void update(const BeginOfRun* run) override;
58  void update(const BeginOfEvent* evt) override;
59  void update(const EndOfEvent* evt) override;
60  void update(const G4Step* step) override;
61 
62  // analysis-related stuff
63  std::vector<int> layerGrouping(int);
64  std::vector<int> towersToAdd(int centre, int nadd);
65  void fill(const EndOfEvent* ev);
66  void qieAnalysis(CLHEP::HepRandomEngine*);
67  void layerAnalysis();
68  double timeOfFlight(int det, int layer, double eta);
69 
70 private:
71  // Qie Analysis
72  std::unique_ptr<HcalQie> myqie_;
73  int addTower_;
74 
75  // Private Tuples
76  std::unique_ptr<HcalTestHistoClass> tuples_;
77 
78  // to read from ParameterSet
80  double eta0_, phi0_;
82  const std::vector<std::string> names_;
83  //Keep parameters to instantiate HcalTestHistoClass later
85 
86  // Numbering scheme
88  std::unique_ptr<HcalNumberingFromDDD> numberingFromDDD_;
90  std::unique_ptr<HcalTestNumberingScheme> org_;
91 
92  // Hits for qie analysis
93  std::vector<CaloHit> caloHitCache_;
94  std::vector<int> group_, tower_;
96 
97  // some private members for ananlysis
98  unsigned int count_;
100  double edepHO_, edepl_[20];
101  double mudist_[20]; // Distance of muon from central part
102 };
103 
105  : addTower_(3),
106  m_Anal(p.getParameter<edm::ParameterSet>("HcalTestAnalysis")),
107  eta0_(m_Anal.getParameter<double>("Eta0")),
108  phi0_(m_Anal.getParameter<double>("Phi0")),
109  laygroup_(m_Anal.getParameter<int>("LayerGrouping")),
110  centralTower_(m_Anal.getParameter<int>("CentralTower")),
111  names_(m_Anal.getParameter<std::vector<std::string> >("Names")),
112  fileName_(m_Anal.getParameter<std::string>("FileName")),
113  hcons_(nullptr) {
114  org_.reset(nullptr);
115 
116  tuples_.reset(nullptr);
117  numberingFromDDD_.reset(nullptr);
118  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: Initialised as observer of begin/end events"
119  << " and of G4step";
120 
121  count_ = 0;
123  nGroup_ = 0;
124  for (unsigned int i = 0; i < group_.size(); i++)
125  if (group_[i] > nGroup_)
126  nGroup_ = group_[i];
128  nTower_ = tower_.size() / 2;
129 
130  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: initialised for " << nGroup_ << " Longitudinal groups and "
131  << nTower_ << " towers";
132 
133  // qie
134  myqie_ = std::make_unique<HcalQie>(p);
135 
136  produces<HcalTestHistoClass>();
137 }
138 
140  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis: --------> Total number of selected entries : " << count_;
141  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis: Pointers:: Numbering Scheme " << org_.get();
142 }
143 
146  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::Initialize ESGetToken for HcalDDDSimConstants";
147 }
148 
150  e.put(std::move(tuples_));
151  tuples_.reset(nullptr);
152 }
153 
154 std::vector<int> HcalTestAnalysis::layerGrouping(int group) {
155  std::vector<int> temp(19);
156  if (group <= 1) {
157  int grp[19] = {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2};
158  for (int i = 0; i < 19; i++)
159  temp[i] = grp[i];
160  } else if (group == 2) {
161  int grp[19] = {1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19};
162  for (int i = 0; i < 19; i++)
163  temp[i] = grp[i];
164  } else if (group == 3) {
165  int grp[19] = {1, 1, 2, 2, 2, 3, 3, 3, 4, 4, 4, 5, 5, 5, 6, 6, 6, 7, 7};
166  for (int i = 0; i < 19; i++)
167  temp[i] = grp[i];
168  } else if (group == 4) {
169  int grp[19] = {1, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 6, 6, 6, 7, 7};
170  for (int i = 0; i < 19; i++)
171  temp[i] = grp[i];
172  } else {
173  int grp[19] = {1, 1, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 6, 6, 6, 6, 6, 7, 7};
174  for (int i = 0; i < 19; i++)
175  temp[i] = grp[i];
176  }
177 
178  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: Layer Grouping ";
179  for (int i = 0; i < 19; i++)
180  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis: Group[" << i << "] = " << temp[i];
181  return temp;
182 }
183 
184 std::vector<int> HcalTestAnalysis::towersToAdd(int centre, int nadd) {
185  int etac = (centre / 100) % 100;
186  int phic = (centre % 100);
187  int etamin, etamax, phimin, phimax;
188  if (etac > 0) {
189  etamin = etac - nadd;
190  etamax = etac + nadd;
191  } else {
192  etamin = etac;
193  etamax = etac;
194  }
195  if (phic > 0) {
196  phimin = phic - nadd;
197  phimax = phic + nadd;
198  } else {
199  phimin = phic;
200  phimax = phic;
201  }
202 
203  int nbuf, kount = 0;
204  nbuf = (etamax - etamin + 1) * (phimax - phimin + 1);
205  std::vector<int> temp(2 * nbuf);
206  for (int eta = etamin; eta <= etamax; eta++) {
207  for (int phi = phimin; phi <= phimax; phi++) {
208  temp[kount] = (eta * 100 + phi);
209  temp[kount + nbuf] = std::max(abs(eta - etac), abs(phi - phic));
210  kount++;
211  }
212  }
213 
214  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: Towers to be considered for Central " << centre << " and " << nadd
215  << " on either side";
216  for (int i = 0; i < nbuf; i++)
217  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis: Tower[" << std::setw(3) << i << "] " << temp[i] << " "
218  << temp[nbuf + i];
219  return temp;
220 }
221 
222 //==================================================================== per JOB
223 
225  // Numbering From DDD
226  hcons_ = &es.getData(ddconsToken_);
227  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: Initialise HcalNumberingFromDDD for " << names_[0];
228  numberingFromDDD_ = std::make_unique<HcalNumberingFromDDD>(hcons_);
229 
230  // Numbering scheme
231  org_ = std::make_unique<HcalTestNumberingScheme>(false);
232 }
233 
234 //==================================================================== per RUN
236  int irun = (*run)()->GetRunID();
237  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: Begin of Run = " << irun;
238 
239  bool loop = true, eta = true, phi = true;
240  int etac = (centralTower_ / 100) % 100;
241  if (etac == 0) {
242  etac = 1;
243  eta = false;
244  }
245  int phic = (centralTower_ % 100);
246  if (phic == 0) {
247  phic = 1;
248  phi = false;
249  }
250  int idet = static_cast<int>(HcalBarrel);
251  while (loop) {
252  HcalCellType::HcalCell tmp = hcons_->cell(idet, 1, 1, etac, phic);
253  if (tmp.ok) {
254  if (eta)
255  eta0_ = tmp.eta;
256  if (phi)
257  phi0_ = tmp.phi;
258  loop = false;
259  } else if (idet == static_cast<int>(HcalBarrel)) {
260  idet = static_cast<int>(HcalEndcap);
261  } else if (idet == static_cast<int>(HcalEndcap)) {
262  idet = static_cast<int>(HcalForward);
263  } else {
264  loop = false;
265  }
266  }
267 
268  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: Central Tower " << centralTower_
269  << " corresponds to eta0 = " << eta0_ << " phi0 = " << phi0_;
270 
271  std::string sdname = names_[0];
272  G4SDManager* sd = G4SDManager::GetSDMpointerIfExist();
273  if (sd != nullptr) {
274  G4VSensitiveDetector* aSD = sd->FindSensitiveDetector(sdname);
275  if (aSD == nullptr) {
276  edm::LogWarning("HcalSim") << "HcalTestAnalysis::beginOfRun: No SD with "
277  << "name " << sdname << " in this Setup";
278  } else {
279  HCalSD* theCaloSD = dynamic_cast<HCalSD*>(aSD);
280  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::beginOfRun: Finds SD with name " << theCaloSD->GetName()
281  << " in this Setup";
282  if (org_.get()) {
283  theCaloSD->setNumberingScheme(org_.get());
284  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::beginOfRun: set a new numbering scheme";
285  }
286  }
287  } else {
288  edm::LogWarning("HcalSim") << "HcalTestAnalysis::beginOfRun: Could not get"
289  << " SD Manager!";
290  }
291 }
292 
293 //=================================================================== per EVENT
295  // create tuple object
296  tuples_ = std::make_unique<HcalTestHistoClass>();
297  // Reset counters
298  tuples_->setCounters();
299 
300  int i = 0;
301  edepEB_ = edepEE_ = edepHB_ = edepHE_ = edepHO_ = 0.;
302  for (i = 0; i < 20; i++)
303  edepl_[i] = 0.;
304  for (i = 0; i < 20; i++)
305  mudist_[i] = -1.;
306 
307  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis: Begin of event = " << (*evt)()->GetEventID();
308 }
309 
310 //=================================================================== each STEP
311 void HcalTestAnalysis::update(const G4Step* aStep) {
312  if (aStep != nullptr) {
313  G4VPhysicalVolume* curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
314  G4String name = curPV->GetName();
315  name.assign(name, 0, 3);
316  double edeposit = aStep->GetTotalEnergyDeposit();
317  int layer = -1;
318  if (name == "EBR") {
319  edepEB_ += edeposit;
320  } else if (name == "EFR") {
321  edepEE_ += edeposit;
322  } else if (name == "HBS") {
323  layer = (curPV->GetCopyNo() / 10) % 100;
324  if (layer >= 0 && layer < 17) {
325  edepHB_ += edeposit;
326  } else {
327  edm::LogWarning("HcalSim") << "HcalTestAnalysis::Error in HB " << curPV->GetName() << curPV->GetCopyNo();
328  layer = -1;
329  }
330  } else if (name == "HES") {
331  layer = (curPV->GetCopyNo() / 10) % 100;
332  if (layer >= 0 && layer < 19) {
333  edepHE_ += edeposit;
334  } else {
335  edm::LogWarning("HcalSim") << "HcalTestAnalysis::Error in HE " << curPV->GetName() << curPV->GetCopyNo();
336  layer = -1;
337  }
338  } else if (name == "HTS") {
339  layer = (curPV->GetCopyNo() / 10) % 100;
340  if (layer >= 17 && layer < 20) {
341  edepHO_ += edeposit;
342  } else {
343  edm::LogWarning("HcalSim") << "HcalTestAnalysis::Error in HO " << curPV->GetName() << curPV->GetCopyNo();
344  layer = -1;
345  }
346  }
347  if (layer >= 0 && layer < 20) {
348  edepl_[layer] += edeposit;
349 
350  // Calculate the distance if it is a muon
351  G4String part = aStep->GetTrack()->GetDefinition()->GetParticleName();
352  if ((part == "mu-" || part == "mu+") && mudist_[layer] < 0) {
353  math::XYZPoint pos(aStep->GetPreStepPoint()->GetPosition().x(),
354  aStep->GetPreStepPoint()->GetPosition().y(),
355  aStep->GetPreStepPoint()->GetPosition().z());
356  double theta = pos.theta();
357  double eta = -log(tan(theta * 0.5));
358  double phi = pos.phi();
359  double dist = sqrt((eta - eta0_) * (eta - eta0_) + (phi - phi0_) * (phi - phi0_));
360  mudist_[layer] = dist * std::sqrt(pos.perp2());
361  }
362  }
363 
364  if (layer >= 0 && layer < 20) {
365  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: G4Step: " << name << " Layer " << std::setw(3) << layer
366  << " Edep " << std::setw(6) << edeposit / MeV << " MeV";
367  }
368  } else {
369  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: G4Step: Null Step";
370  }
371 }
372 
373 //================================================================ End of EVENT
375  ++count_;
376  // Fill event input
377  fill(evt);
378  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: --- after Fill";
379 
380  // Qie analysis
381  CLHEP::HepRandomEngine* engine = G4Random::getTheEngine();
382  qieAnalysis(engine);
383  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: --- after QieAnalysis";
384 
385  // Layers tuples filling
386  layerAnalysis();
387  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis:: --- after LayerAnalysis";
388 }
389 
390 //---------------------------------------------------
392  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis: Fill event " << (*evt)()->GetEventID();
393 
394  // access to the G4 hit collections
395  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
396 
397  int nhc = 0, neb = 0, nef = 0, j = 0;
398  caloHitCache_.erase(caloHitCache_.begin(), caloHitCache_.end());
399 
400  // Hcal
401  int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names_[0]);
402  CaloG4HitCollection* theHCHC = (CaloG4HitCollection*)allHC->GetHC(HCHCid);
403  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis :: Hit Collection for " << names_[0] << " of ID " << HCHCid
404  << " is obtained at " << theHCHC;
405  int hchc_entries = theHCHC->entries();
406  if (HCHCid >= 0 && theHCHC != nullptr) {
407  for (j = 0; j < hchc_entries; j++) {
408  CaloG4Hit* aHit = (*theHCHC)[j];
409 
410  double e = aHit->getEnergyDeposit() / GeV;
411  double time = aHit->getTimeSlice();
412 
413  math::XYZPoint pos = aHit->getPosition();
414  double theta = pos.theta();
415  double eta = -log(tan(theta * 0.5));
416  double phi = pos.phi();
417 
418  uint32_t unitID = aHit->getUnitID();
419  int subdet, zside, layer, etaIndex, phiIndex, lay;
420  org_->unpackHcalIndex(unitID, subdet, zside, layer, etaIndex, phiIndex, lay);
421  double jitter = time - timeOfFlight(subdet, lay, eta);
422  if (jitter < 0)
423  jitter = 0;
424  CaloHit hit(subdet, lay, e, eta, phi, jitter, unitID);
425  caloHitCache_.push_back(hit);
426  nhc++;
427 
428  std::string det = "HB";
429  if (subdet == static_cast<int>(HcalForward)) {
430  det = "HF";
431  } else if (subdet == static_cast<int>(HcalEndcap)) {
432  if (etaIndex <= 20) {
433  det = "HES";
434  } else {
435  det = "HED";
436  }
437  }
438 
439  edm::LogVerbatim("HcalSim") << "HcalTest: " << det << " layer " << std::setw(2) << layer << " time "
440  << std::setw(6) << time << " theta " << std::setw(8) << theta << " eta "
441  << std::setw(8) << eta << " phi " << std::setw(8) << phi << " e " << std::setw(8)
442  << e;
443  }
444  }
445  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::HCAL hits : " << nhc;
446 
447  // EB
448  int EBHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names_[1]);
449  CaloG4HitCollection* theEBHC = (CaloG4HitCollection*)allHC->GetHC(EBHCid);
450  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis :: Hit Collection for " << names_[1] << " of ID " << EBHCid
451  << " is obtained at " << theEBHC;
452  int ebhc_entries = theEBHC->entries();
453  if (EBHCid >= 0 && theEBHC != nullptr) {
454  for (j = 0; j < ebhc_entries; j++) {
455  CaloG4Hit* aHit = (*theEBHC)[j];
456 
457  double e = aHit->getEnergyDeposit() / GeV;
458  double time = aHit->getTimeSlice();
459  std::string det = "EB";
460 
461  math::XYZPoint pos = aHit->getPosition();
462  double theta = pos.theta();
463  double eta = -log(tan(theta / 2.));
464  double phi = pos.phi();
465 
467  uint32_t unitID = org_->getUnitID(id);
468  int subdet, zside, layer, ieta, iphi, lay;
469  org_->unpackHcalIndex(unitID, subdet, zside, layer, ieta, iphi, lay);
470  subdet = 10;
471  layer = 0;
472  unitID = org_->packHcalIndex(subdet, zside, layer, ieta, iphi, lay);
473  CaloHit hit(subdet, lay, e, eta, phi, time, unitID);
474  caloHitCache_.push_back(hit);
475  neb++;
476  edm::LogVerbatim("HcalSim") << "HcalTest: " << det << " layer " << std::setw(2) << layer << " time "
477  << std::setw(6) << time << " theta " << std::setw(8) << theta << " eta "
478  << std::setw(8) << eta << " phi " << std::setw(8) << phi << " e " << std::setw(8)
479  << e;
480  }
481  }
482  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::EB hits : " << neb;
483 
484  // EE
485  int EEHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names_[2]);
486  CaloG4HitCollection* theEEHC = (CaloG4HitCollection*)allHC->GetHC(EEHCid);
487  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis :: Hit Collection for " << names_[2] << " of ID " << EEHCid
488  << " is obtained at " << theEEHC;
489  int eehc_entries = theEEHC->entries();
490  if (EEHCid >= 0 && theEEHC != nullptr) {
491  for (j = 0; j < eehc_entries; j++) {
492  CaloG4Hit* aHit = (*theEEHC)[j];
493 
494  double e = aHit->getEnergyDeposit() / GeV;
495  double time = aHit->getTimeSlice();
496  std::string det = "EE";
497 
498  math::XYZPoint pos = aHit->getPosition();
499  double theta = pos.theta();
500  double eta = -log(tan(theta / 2.));
501  double phi = pos.phi();
502 
504  uint32_t unitID = org_->getUnitID(id);
505  int subdet, zside, layer, ieta, iphi, lay;
506  org_->unpackHcalIndex(unitID, subdet, zside, layer, ieta, iphi, lay);
507  subdet = 11;
508  layer = 0;
509  unitID = org_->packHcalIndex(subdet, zside, layer, ieta, iphi, lay);
510  CaloHit hit(subdet, lay, e, eta, phi, time, unitID);
511  caloHitCache_.push_back(hit);
512  nef++;
513  edm::LogVerbatim("HcalSim") << "HcalTest: " << det << " layer " << std::setw(2) << layer << " time "
514  << std::setw(6) << time << " theta " << std::setw(8) << theta << " eta "
515  << std::setw(8) << eta << " phi " << std::setw(8) << phi << " e " << std::setw(8)
516  << e;
517  }
518  }
519  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::EE hits : " << nef;
520 }
521 
522 //-----------------------------------------------------------------------------
523 void HcalTestAnalysis::qieAnalysis(CLHEP::HepRandomEngine* engine) {
524  //Fill tuple with hit information
525  int hittot = caloHitCache_.size();
526  tuples_->fillHits(caloHitCache_);
527 
528  //Get the index of the central tower
530  uint32_t unitID = org_->getUnitID(id);
531  int subdet, zside, layer, ieta, iphi, lay;
532  org_->unpackHcalIndex(unitID, subdet, zside, layer, ieta, iphi, lay);
533  int laymax = 0;
534  std::string det = "Unknown";
535  if (subdet == static_cast<int>(HcalBarrel)) {
536  laymax = 4;
537  det = "HB";
538  } else if (subdet == static_cast<int>(HcalEndcap)) {
539  laymax = 2;
540  det = "HES";
541  }
542  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::Qie: " << det << " Eta " << ieta << " Phi " << iphi << " Laymax "
543  << laymax << " Hits " << hittot;
544 
545  if (laymax > 0 && hittot > 0) {
546  std::vector<CaloHit> hits(hittot);
547  std::vector<double> eqielay(80, 0.0), esimlay(80, 0.0), esimtot(4, 0.0);
548  std::vector<double> eqietow(200, 0.0), esimtow(200, 0.0), eqietot(4, 0.0);
549  int etac = (centralTower_ / 100) % 100;
550  int phic = (centralTower_ % 100);
551 
552  for (int layr = 0; layr < nGroup_; layr++) {
553  /*
554  int layx, layy=20;
555  for (int i=0; i<20; i++)
556  if (group_[i] == layr+1 && i < layy) layy = i+1;
557  if (subdet == static_cast<int>(HcalBarrel)) {
558  if (layy < 2) layx = 0;
559  else if (layy < 17) layx = 1;
560  else if (layy == 17) layx = 2;
561  else layx = 3;
562  } else {
563  if (layy < 2) layx = 0;
564  else layx = 1;
565  }
566  */
567  for (int it = 0; it < nTower_; it++) {
568  int nhit = 0;
569  double esim = 0;
570  for (int k1 = 0; k1 < hittot; k1++) {
571  CaloHit hit = caloHitCache_[k1];
572  int subdetc = hit.det();
573  int layer = hit.layer();
574  int group = 0;
575  if (layer > 0 && layer < 20)
576  group = group_[layer];
577  if (subdetc == subdet && group == layr + 1) {
578  int zsidec, ietac, iphic, idx;
579  unitID = hit.id();
580  org_->unpackHcalIndex(unitID, subdetc, zsidec, layer, ietac, iphic, lay);
581  if (etac > 0 && phic > 0) {
582  idx = ietac * 100 + iphic;
583  } else if (etac > 0) {
584  idx = ietac * 100;
585  } else if (phic > 0) {
586  idx = iphic;
587  } else {
588  idx = 0;
589  }
590  if (zsidec == zside && idx == tower_[it]) {
591  hits[nhit] = hit;
592  edm::LogVerbatim("HcalSim") << "HcalTest: Hit " << nhit << " " << hit;
593  nhit++;
594  esim += hit.e();
595  }
596  }
597  }
598 
599  std::vector<int> cd = myqie_->getCode(nhit, hits, engine);
600  double eqie = myqie_->getEnergy(cd);
601 
602  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::Qie: Energy in layer " << layr << " Sim " << esim
603  << " After QIE " << eqie;
604  for (int i = 0; i < 4; i++) {
605  if (tower_[nTower_ + it] <= i) {
606  esimtot[i] += esim;
607  eqietot[i] += eqie;
608  esimlay[20 * i + layr] += esim;
609  eqielay[20 * i + layr] += eqie;
610  esimtow[50 * i + it] += esim;
611  eqietow[50 * i + it] += eqie;
612  }
613  }
614  }
615  }
616  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::Qie: Total energy " << esimtot[3] << " (SimHit) " << eqietot[3]
617  << " (After QIE)";
618 
619  std::vector<double> latphi(10);
620  int nt = 2 * addTower_ + 1;
621  for (int it = 0; it < nt; it++)
622  latphi[it] = it - addTower_;
623  for (int i = 0; i < 4; i++) {
624  double scals = 1, scalq = 1;
625  std::vector<double> latfs(10, 0.), latfq(10, 0.), longs(20), longq(20);
626  if (esimtot[i] > 0)
627  scals = 1. / esimtot[i];
628  if (eqietot[i] > 0)
629  scalq = 1. / eqietot[i];
630  for (int it = 0; it < nTower_; it++) {
631  int phib = it % nt;
632  latfs[phib] += scals * esimtow[50 * i + it];
633  latfq[phib] += scalq * eqietow[50 * i + it];
634  }
635  for (int layr = 0; layr <= nGroup_; layr++) {
636  longs[layr] = scals * esimlay[20 * i + layr];
637  longq[layr] = scalq * eqielay[20 * i + layr];
638  }
639  tuples_->fillQie(i, esimtot[i], eqietot[i], nGroup_, longs, longq, nt, latphi, latfs, latfq);
640  }
641  }
642 }
643 
644 //---------------------------------------------------
646  int i = 0;
647  edm::LogVerbatim("HcalSim") << "\n ===>>> HcalTestAnalysis: Energy deposit "
648  << "\n at EB : " << std::setw(6) << edepEB_ / MeV << "\n at EE : " << std::setw(6)
649  << edepEE_ / MeV << "\n at HB : " << std::setw(6) << edepHB_ / MeV
650  << "\n at HE : " << std::setw(6) << edepHE_ / MeV << "\n at HO : " << std::setw(6)
651  << edepHO_ / MeV << "\n ---- HcalTestAnalysis: Energy deposit in Layers";
652  for (i = 0; i < 20; i++)
653  edm::LogVerbatim("HcalSim") << " Layer " << std::setw(2) << i << " E " << std::setw(8) << edepl_[i] / MeV << " MeV";
654 
655  tuples_->fillLayers(edepl_, edepHO_, edepHB_ + edepHE_, mudist_);
656 }
657 
658 //---------------------------------------------------
659 double HcalTestAnalysis::timeOfFlight(int det, int layer, double eta) {
660  double theta = 2.0 * atan(exp(-eta));
661  double dist = 0.;
662  if (det == static_cast<int>(HcalBarrel)) {
663  const double rLay[19] = {1836.0,
664  1902.0,
665  1962.0,
666  2022.0,
667  2082.0,
668  2142.0,
669  2202.0,
670  2262.0,
671  2322.0,
672  2382.0,
673  2448.0,
674  2514.0,
675  2580.0,
676  2646.0,
677  2712.0,
678  2776.0,
679  2862.5,
680  3847.0,
681  4052.0};
682  if (layer > 0 && layer < 20)
683  dist += rLay[layer - 1] * mm / sin(theta);
684  } else {
685  const double zLay[19] = {4034.0,
686  4032.0,
687  4123.0,
688  4210.0,
689  4297.0,
690  4384.0,
691  4471.0,
692  4558.0,
693  4645.0,
694  4732.0,
695  4819.0,
696  4906.0,
697  4993.0,
698  5080.0,
699  5167.0,
700  5254.0,
701  5341.0,
702  5428.0,
703  5515.0};
704  if (layer > 0 && layer < 20)
705  dist += zLay[layer - 1] * mm / cos(theta);
706  }
707  double tmp = dist / c_light / ns;
708  edm::LogVerbatim("HcalSim") << "HcalTestAnalysis::timeOfFlight " << tmp << " for det/lay " << det << " " << layer
709  << " eta/theta " << eta << " " << theta / deg << " dist " << dist;
710  return tmp;
711 }
712 
715 
HcalCellType::HcalCell cell(const int &det, const int &zside, const int &depth, const int &etaR, const int &iphi) const
Log< level::Info, true > LogVerbatim
#define DEFINE_SIMWATCHER(type)
const std::vector< std::string > names_
std::unique_ptr< HcalTestHistoClass > tuples_
std::vector< int > layerGrouping(int)
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const edm::ParameterSet m_Anal
int zside(DetId const &)
void beginRun(edm::EventSetup const &) override
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:544
std::vector< int > towersToAdd(int centre, int nadd)
constexpr std::array< uint8_t, layerIndexSize > layer
void produce(edm::Event &, const edm::EventSetup &) override
double timeOfFlight(int det, int layer, double eta)
Definition: HCalSD.h:38
math::XYZPoint getPosition() const
Definition: CaloG4Hit.h:52
const std::string fileName_
void fill(const EndOfEvent *ev)
std::vector< int > tower_
HcalTestAnalysis(const edm::ParameterSet &p)
T sqrt(T t)
Definition: SSEVec.h:19
const HcalDDDSimConstants * hcons_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
unsigned int count_
~HcalTestAnalysis() override
std::unique_ptr< HcalQie > myqie_
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getData(T &iHolder) const
Definition: EventSetup.h:122
int nt
Definition: AMPTWrapper.h:42
unsigned int id
void registerConsumes(edm::ConsumesCollector) override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::unique_ptr< HcalTestNumberingScheme > org_
part
Definition: HCALResponse.h:20
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
std::vector< CaloHit > caloHitCache_
HLT enums.
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
double getTimeSlice() const
Definition: CaloG4Hit.h:67
edm::ESGetToken< HcalDDDSimConstants, HcalSimNumberingRecord > ddconsToken_
step
Definition: StallMonitor.cc:98
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
void qieAnalysis(CLHEP::HepRandomEngine *)
Log< level::Warning, false > LogWarning
tmp
align.sh
Definition: createJobs.py:716
Geom::Theta< T > theta() const
std::vector< int > group_
def move(src, dest)
Definition: eostools.py:511