CMS 3D CMS Logo

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