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