CMS 3D CMS Logo

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