CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalTB02Analysis.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestBeam
4 // Class : HcalTB02Analysis
5 //
6 // Implementation:
7 // Main analysis class for Hcal Test Beam 2002 Analysis
8 //
9 // Original Author:
10 // Created: Sun May 21 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 #include <cmath>
15 #include <iostream>
16 #include <iomanip>
17 
18 // user include files
20 
21 // to retreive hits
27 
32 #include "CLHEP/Random/RandGaussQ.h"
33 
34 #include "G4SDManager.hh"
35 #include "G4VProcess.hh"
36 #include "G4HCofThisEvent.hh"
37 #include "CLHEP/Units/GlobalSystemOfUnits.h"
38 #include "CLHEP/Units/GlobalPhysicalConstants.h"
39 #include "CLHEP/Random/Random.h"
40 
41 namespace CLHEP {
42  class HepRandomEngine;
43 }
44 
45 //
46 // constructors and destructor
47 //
48 
50 
51  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTB02Analysis");
52  hcalOnly = m_Anal.getUntrackedParameter<bool>("HcalClusterOnly",true);
53  names = m_Anal.getParameter<std::vector<std::string> >("Names");
54 
55  produces<HcalTB02HistoClass>();
56 
57  edm::LogInfo("HcalTBSim") << "HcalTB02Analysis:: Initialised as observer of "
58  << "BeginOfJob/BeginOfEvent/EndOfEvent with "
59  << "Parameter values:\n \thcalOnly = " << hcalOnly;
60 
61  histo = new HcalTB02Histo(m_Anal);
62 }
63 
65 
66  finish();
67 
68  if (histo) {
69  delete histo;
70  histo = 0;
71  }
72  edm::LogInfo("HcalTBSim") << "HcalTB02Analysis is deleting";
73 }
74 
75 //
76 // member functions
77 //
78 
80 
81  std::auto_ptr<HcalTB02HistoClass> product(new HcalTB02HistoClass);
82  fillEvent(*product);
83  e.put(product);
84 }
85 
87 
88  edm::LogInfo("HcalTBSim") << "HcalTB02Analysis: =====> Begin of event = "
89  << (*evt) ()->GetEventID();
90  clear();
91 }
92 
94 
95  CLHEP::HepRandomEngine* engine = CLHEP::HepRandom::getTheEngine();
96  CLHEP::RandGaussQ randGauss(*engine);
97 
98  // Look for the Hit Collection
99  LogDebug("HcalTBSim") << "HcalTB02Analysis::Fill event "
100  << (*evt)()->GetEventID();
101 
102  // access to the G4 hit collections
103  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
104  int ihit = 0;
105 
106  // HCAL
107  std::string sd = names[0];
108  int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(sd);
109  CaloG4HitCollection* theHCHC = (CaloG4HitCollection*) allHC->GetHC(HCHCid);
111  LogDebug("HcalTBSim") << "HcalTB02Analysis :: Hit Collection for " << sd
112  << " of ID " << HCHCid << " is obtained at " <<theHCHC;
113 
114  int nentries = 0;
115  nentries = theHCHC->entries();
116  if (nentries==0) return;
117 
118  int xentries = 0;
119  int XTALSid=0;
120  CaloG4HitCollection* theXTHC=0;
121 
122  if (!hcalOnly) {
123  // XTALS
124  sd = names[1];
125  XTALSid = G4SDManager::GetSDMpointer()->GetCollectionID(sd);
126  // assert (XTALSid != 0);
127  theXTHC = (CaloG4HitCollection*) allHC->GetHC(XTALSid);
128  // assert (theXTHC != 0);
129  //HcalTB02XtalNumberingScheme *xorg = new HcalTB02XtalNumberingScheme();
130  LogDebug("HcalTBSim") << "HcalTB02Analysis :: Hit Collection for " << sd
131  << " of ID " << XTALSid << " is obtained at "
132  << theXTHC;
133  xentries = theXTHC->entries();
134  if (xentries==0) return;
135  }
136 
137  LogDebug("HcalTBSim") << "HcalTB02Analysis :: There are " << nentries
138  << " HCal hits, and" << xentries << " xtal hits";
139 
140  float ETot=0., xETot=0.;
141  //float maxE = 0.;
142  //int maxI=0,
143  int scintID=0, xtalID=0;
144 
145  // HCAL
146 
147  if (HCHCid >= 0 && theHCHC > 0) {
148  for ( ihit = 0; ihit < nentries; ihit++) {
149 
150  CaloG4Hit* aHit = (*theHCHC)[ihit];
151  scintID = aHit->getUnitID();
152  int layer = org->getlayerID(scintID);
153  float enEm = aHit->getEM();
154  float enhad = aHit->getHadr();
155 
156  if (layer==0) {
157  enEm =enEm/2.;
158  enhad=enhad/2.;
159  }
160 
161  energyInScints[scintID]+= enEm + enhad;
162  primaries[aHit->getTrackID()]+= enEm + enhad;
163  float time = aHit->getTimeSliceID();
164  if (time > maxTime) maxTime=(int)time;
165  histo->fillAllTime(time);
166 
167  }
168 
169  //
170  // Profile
171  //
172 
173  float TowerEne[8][18], TowerEneCF[8][18], LayerEne[19], EnRing[100];
174  for (int iphi=0 ; iphi<8; iphi++) {
175  for (int jeta=0 ; jeta<18; jeta++) {
176  TowerEne[iphi][jeta]=0.;
177  TowerEneCF[iphi][jeta]=0.;
178  }
179  }
180 
181  for (int ilayer=0; ilayer<19; ilayer++) LayerEne[ilayer]=0.;
182  for (int iring=0; iring<100; iring++) EnRing[iring]=0.;
183 
184  for (std::map<int,float>::iterator is = energyInScints.begin();
185  is!= energyInScints.end(); is++) {
186 
187  ETot = (*is).second;
188 
189  int layer = org->getlayerID((*is).first);
190 
191  if ( (layer!=17) && (layer!=18) ) {
192 
193  float eta = org->getetaID((*is).first);
194  float phi = org->getphiID((*is).first);
195 
196  SEnergy += ETot;
197  TowerEne[(int)phi][(int)eta] += ETot;
198 
199  TowerEneCF[(int)phi][(int)eta] += ETot*(1.+ 0.1*randGauss.fire() );
200  double dR=0.08727*std::sqrt( (eta-8.)*(eta-8.) + (phi-3.)*(phi-3.) );
201  EnRing[(int)(dR/0.01)] += ETot;
202  }
203 
204  LayerEne[layer] += ETot;
205 
206  }
207 
208  for (int ilayer=0 ; ilayer<19 ; ilayer++) {
209  histo->fillProfile(ilayer,LayerEne[ilayer]/GeV);
210  }
211 
212  for (int iring=0; iring<100; iring++)
213  histo->fillTransProf(0.01*iring+0.005,EnRing[iring]/GeV);
214 
215  for (int iphi=0 ; iphi<8; iphi++) {
216  for (int jeta=0 ; jeta<18; jeta++) {
217 
218  //SEnergyN += TowerEneCF[iphi][jeta] + 3.2*randGauss.fire(); // LHEP
219  SEnergyN += TowerEneCF[iphi][jeta] + 3.*randGauss.fire(); // QGSP
220 
221  //double dR=0.08727*sqrt( (jeta-8.)*(jeta-8.)+(iphi-3.)*(iphi-3.) );
222  //cout.testOut << " phi= " << iphi << " eta= " << jeta
223  // << " TowerEne[iphi,jeta]= " << TowerEne[iphi][jeta]
224  // << "dR= " << dR << endl;
225 
226  //double Rand = 3.2*randGauss.fire(); // LHEP
227  double Rand = 3.*randGauss.fire(); // QGSP
228 
229  if ( (iphi>=0) && (iphi<7) ) {
230  if ( (jeta>=5) && (jeta<12) ) {
231 
232  E7x7Matrix += TowerEne[iphi][jeta];
233  E7x7MatrixN += TowerEneCF[iphi][jeta] + Rand;
234 
235  if ( (iphi>=1) && (iphi<6) ) {
236  if ( (jeta>=6) && (jeta<11) ) {
237 
238  E5x5Matrix += TowerEne[iphi][jeta];
239  E5x5MatrixN += TowerEneCF[iphi][jeta] + Rand;
240 
241  }
242  }
243 
244  }
245  }
246 
247  }
248  }
249 
250  //
251  // Find Primary info:
252  //
253  int trackID = 0;
254  G4PrimaryParticle* thePrim=0;
255  G4int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
256  LogDebug("HcalTBSim") << "HcalTB02Analysis :: Event has " << nvertex
257  << " vertex";
258  if (nvertex==0)
259  edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
260  << "ERROR: no vertex";
261 
262  for (int i = 0 ; i<nvertex; i++) {
263 
264  G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
265  if (avertex == 0) {
266  edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
267  << "ERROR: pointer to vertex = 0";
268  } else {
269  int npart = avertex->GetNumberOfParticle();
270  LogDebug("HcalTBSim") << "HcalTB02Analysis::Vertex number :" << i
271  << " with " << npart << " particles";
272  if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
273  }
274  }
275 
276  double px=0.,py=0.,pz=0.;
277 
278  if (thePrim != 0) {
279  px = thePrim->GetPx();
280  py = thePrim->GetPy();
281  pz = thePrim->GetPz();
282  pInit = std::sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
283  if (pInit==0) {
284  edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
285  << " ERROR: primary has p=0 ";
286  } else {
287  float costheta = pz/pInit;
288  float theta = acos(std::min(std::max(costheta,float(-1.)),float(1.)));
289  eta = -log(tan(theta/2));
290  if (px != 0) phi = atan(py/px);
291  }
292  particleType = thePrim->GetPDGcode();
293  } else {
294  LogDebug("HcalTBSim") << "HcalTB02Analysis:: End Of Event ERROR: could"
295  << " not find primary ";
296  }
297 
298  CaloG4Hit* firstHit =(*theHCHC)[0];
299  incidentEnergy = (firstHit->getIncidentEnergy()/GeV);
300 
301  }// number of Hits > 0
302 
303  if (!hcalOnly) {
304 
305  // XTALS
306 
307  if (XTALSid >= 0 && theXTHC > 0) {
308  for (int xihit = 0; xihit < xentries; xihit++) {
309 
310  CaloG4Hit* xaHit = (*theXTHC)[xihit];
311 
312  float xenEm = xaHit->getEM();
313  float xenhad = xaHit->getHadr();
314  xtalID = xaHit->getUnitID();
315 
316  energyInCrystals[xtalID]+= xenEm + xenhad;
317  }
318 
319  float xCrysEne[7][7];
320  for (int irow=0 ; irow<7; irow++) {
321  for (int jcol=0 ; jcol<7; jcol++) {
322  xCrysEne[irow][jcol]=0.;
323  }
324  }
325 
326  for (std::map<int,float>::iterator is = energyInCrystals.begin();
327  is!= energyInCrystals.end(); is++) {
328  int xtalID = (*is).first;
329  xETot = (*is).second;
330 
331  int irow = (int)(xtalID/100.);
332  int jcol = (int)(xtalID-100.*irow);
333 
334  xSEnergy += xETot;
335  xCrysEne[irow][jcol] = xETot;
336 
337  float dR=std::sqrt( 0.01619*0.01619*(jcol-3)*(jcol-3) +
338  0.01606*0.01606*(irow-3)*(irow-3) );
339  histo->fillTransProf(dR,xETot*1.05);
340 
341  if ( (irow>0) && (irow<6) ) {
342  if ( (jcol>0) && (jcol<6) ) {
343 
344  xE5x5Matrix += xCrysEne[irow][jcol];
345  xE5x5MatrixN += xCrysEne[irow][jcol] + 108.5*randGauss.fire();
346 
347  if ( (irow>1) && (irow<5) ) {
348  if ( (jcol>1) && (jcol<5) ) {
349  xE3x3Matrix += xCrysEne[irow][jcol];
350  xE3x3MatrixN += xCrysEne[irow][jcol] +108.5*randGauss.fire();
351  }
352  }
353  }
354  }
355 
356  }
357 
358  if (!hcalOnly) {
359  // assert(theXTHC);
360  if ( theXTHC != 0 ) {
361  CaloG4Hit* xfirstHit =(*theXTHC)[0];
362  xIncidentEnergy = xfirstHit->getIncidentEnergy()/GeV;
363  }
364  }
365 
366  }// number of Hits > 0
367 
368  }
369 
370  int iEvt = (*evt)()->GetEventID();
371  if (iEvt < 10)
372  std::cout << " Event " << iEvt << std::endl;
373  else if ((iEvt < 100) && (iEvt%10 == 0))
374  std::cout << " Event " << iEvt << std::endl;
375  else if ((iEvt < 1000) && (iEvt%100 == 0))
376  std::cout << " Event " << iEvt << std::endl;
377  else if ((iEvt < 10000) && (iEvt%1000 == 0))
378  std::cout << " Event " << iEvt << std::endl;
379 }
380 
382 
383  //Beam information
384  product.set_Nprim(float(primaries.size()));
385  product.set_partType(particleType);
386  product.set_Einit(pInit/GeV);
387  product.set_eta(eta);
388  product.set_phi(phi);
389  product.set_Eentry(incidentEnergy);
390 
391  //Calorimeter energy
392  product.set_ETot(SEnergy/GeV );
393  product.set_E7x7(E7x7Matrix/GeV );
394  product.set_E5x5(E5x5Matrix/GeV );
395  product.set_ETotN(SEnergyN/GeV);
396  product.set_E7x7N(E7x7MatrixN/GeV );
397  product.set_E5x5N(E5x5MatrixN/GeV );
398  product.set_NUnit(float(energyInScints.size()));
399  product.set_Ntimesli(float(maxTime));
400 
401  //crystal information
402  product.set_xEentry(xIncidentEnergy);
403  product.set_xNUnit(float(energyInCrystals.size()));
404  product.set_xETot(xSEnergy/GeV);
405  product.set_xETotN(xSEnergyN/GeV);
406  product.set_xE5x5(xE5x5Matrix/GeV);
407  product.set_xE3x3(xE3x3Matrix/GeV);
408  product.set_xE5x5N(xE5x5MatrixN/GeV);
409  product.set_xE3x3N(xE3x3MatrixN/GeV);
410 }
411 
413 
414  primaries.clear();
415  particleType = 0;
416  pInit = eta = phi = incidentEnergy = 0;
417 
419  E7x7MatrixN = E5x5MatrixN = 0;
420  energyInScints.clear();
421  maxTime = 0;
422 
423  xIncidentEnergy = 0;
424  energyInCrystals.clear();
427 }
428 
430 
431  /*
432  //Profile
433  std::ofstream oFile;
434  oFile.open("profile.dat");
435  float st[19] = {0.8,
436  0.4, 0.4, 0.4, 0.4, 0.4,
437  0.4, 0.4, 0.4, 0.4, 0.4,
438  0.4, 0.4, 0.4, 0.4, 0.4,
439  0.8, 1.0, 1.0};
440 
441  //cm of material (brass) in front of scintillator layer i:
442 
443  float w[19] = {7.45, //iron !
444  6.00, 6.00, 6.00, 6.00, 6.00, //brass
445  6.00, 6.00, 6.00, 6.00, 6.60, //brass
446  6.60, 6.60, 6.60, 6.60, 6.60, //brass
447  8.90, 20.65, 19.5}; //brass,iron !
448 
449  for (int ilayer = 0; ilayer<19; ilayer++) {
450 
451  // Histogram mean and sigma calculated from the ROOT histos
452  edm::LogInfo("HcalTBSim") << "Layer number: " << ilayer << " Mean = "
453  << histo->getMean(ilayer) << " sigma = "
454  << histo->getRMS(ilayer) << " LThick= "
455  << w[ilayer] << " SThick= " << st[ilayer];
456 
457  oFile << ilayer << " " << histo->getMean(ilayer) << " "
458  << histo->getRMS(ilayer) << " " << w[ilayer] << " " << st[ilayer]
459  << std::endl;
460 
461  }
462  oFile.close();
463  */
464 }
465 
#define LogDebug(id)
void set_xE3x3N(float v)
void set_E5x5N(float v)
T getParameter(std::string const &) const
void set_Ntimesli(float v)
T getUntrackedParameter(std::string const &, T const &) const
void set_NUnit(float v)
int i
Definition: DBlmapReader.cc:9
#define DEFINE_SIMWATCHER(type)
void set_Eentry(float v)
void set_xEentry(float v)
void set_xETotN(float v)
void set_Nprim(float v)
void set_E7x7(float v)
void set_ETotN(float v)
void fillTransProf(float u, float v)
Geom::Theta< T > theta() const
double getIncidentEnergy() const
Definition: CaloG4Hit.h:65
double npart
Definition: HydjetWrapper.h:44
void set_Einit(float v)
void set_xETot(float v)
void set_E7x7N(float v)
void fillEvent(HcalTB02HistoClass &)
virtual ~HcalTB02Analysis()
std::map< int, float > energyInCrystals
void set_xE3x3(float v)
HcalTB02Analysis(const edm::ParameterSet &p)
HcalTB02Histo * histo
void set_xE5x5N(float v)
const T & max(const T &a, const T &b)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
void fillAllTime(float v)
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
std::map< int, float > primaries
std::map< int, float > energyInScints
int getTrackID() const
Definition: CaloG4Hit.h:68
void update(const BeginOfEvent *evt)
This routine will be called when the appropriate signal arrives.
double sd
void set_E5x5(float v)
virtual void produce(edm::Event &, const edm::EventSetup &)
void set_partType(float v)
int getTimeSliceID() const
Definition: CaloG4Hit.h:71
double getEM() const
Definition: CaloG4Hit.h:59
std::vector< std::string > names
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
void set_ETot(float v)
void fillProfile(int ilayer, float value)
tuple cout
Definition: gather_cfg.py:121
void set_xE5x5(float v)
uint32_t getUnitID() const
Definition: CaloG4Hit.h:69
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double getHadr() const
Definition: CaloG4Hit.h:62
void set_xNUnit(float v)