CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZdcTestAnalysis.cc
Go to the documentation of this file.
1 
5 
8 
9 #include "TFile.h"
10 #include <cmath>
11 #include <iostream>
12 #include <iomanip>
13 
18 };
19 
23 };
24 
26  //constructor
27  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("ZdcTestAnalysis");
28  verbosity = m_Anal.getParameter<int>("Verbosity");
29  doNTzdcstep = m_Anal.getParameter<int>("StepNtupleFlag");
30  doNTzdcevent = m_Anal.getParameter<int>("EventNtupleFlag");
31  stepNtFileName = m_Anal.getParameter<std::string>("StepNtupleFileName");
32  eventNtFileName = m_Anal.getParameter<std::string>("EventNtupleFileName");
33 
34  if (verbosity > 0)
35  std::cout<<std::endl;
36  std::cout<<"============================================================================"<<std::endl;
37  std::cout << "ZdcTestAnalysis:: Initialized as observer"<< std::endl;
38  if (doNTzdcstep > 0){
39  std::cout <<" Step Ntuple will be created"<< std::endl;
40  std::cout <<" Step Ntuple file: "<<stepNtFileName<<std::endl;
41  }
42  if (doNTzdcevent > 0){
43  std::cout <<" Event Ntuple will be created"<< std::endl;
44  std::cout <<" Step Ntuple file: "<<stepNtFileName<<std::endl;
45  }
46  std::cout<<"============================================================================"<<std::endl;
47  std::cout<<std::endl;
48 
49  if (doNTzdcstep > 0)
50  zdcstepntuple =
51  new TNtuple("NTzdcstep","NTzdcstep",
52  "evt:trackid:charge:pdgcode:x:y:z:stepl:stepe:eta:phi:vpx:vpy:vpz:idx:idl:pvtype:ncherphot");
53 
54  if (doNTzdcevent > 0)
56  new TNtuple("NTzdcevent","NTzdcevent",
57  "evt:ihit:fiberid:zside:subdet:layer:fiber:channel:enem:enhad:hitenergy:x:y:z:time:etot");
58 
59  //theZdcSD = new ZdcSD("ZDCHITSB", new ZdcNumberingScheme());
60 }
61 
63  // destructor
64  finish();
65  if (verbosity > 0) {
66  std::cout << std::endl << "ZdcTestAnalysis Dextructor --------> End of ZdcTestAnalysis : "
67  << std::endl << std::endl;
68  }
69 
70  //if (doNTzdcstep > 0)delete zdcstepntuple;
71  //if (doNTzdcevent > 0)delete zdceventntuple;
72 
73  std::cout<<"ZdcTestAnalysis: End of process"<<std::endl;
74 }
75 
76 
78  //job
79  std::cout<<"beggining of job"<<std::endl;;
80 }
81 
82 
83 //==================================================================== per RUN
85  //run
86 
87  std::cout << std::endl << "ZdcTestAnalysis: Begining of Run"<< std::endl;
88  if (doNTzdcstep) {
89  std::cout << "ZDCTestAnalysis: output step file created"<< std::endl;
90  TString stepfilename = stepNtFileName;
91  zdcOutputStepFile = new TFile(stepfilename,"RECREATE");
92 
93  }
94 
95  if (doNTzdcevent) {
96  std::cout << "ZDCTestAnalysis: output event file created"<< std::endl;
97  TString stepfilename = eventNtFileName;
98  zdcOutputEventFile = new TFile(stepfilename,"RECREATE");
99  }
100 
101  eventIndex = 0;
102 }
103 
104 
105 
106 
108  //event
109  std::cout << "ZdcTest: Processing Event Number: "<<eventIndex<< std::endl;
110  eventIndex++;
111  stepIndex = 0;
112 }
113 
114 
115 void ZdcTestAnalysis::update(const G4Step * aStep) {
116  //step;
117  stepIndex++;
118 
119  if (doNTzdcstep) {
120  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
121  // G4StepPoint * postStepPoint= aStep->GetPostStepPoint();
122  G4double stepL = aStep->GetStepLength();
123  G4double stepE = aStep->GetTotalEnergyDeposit();
124 
125  if (verbosity >= 2)
126  std::cout << "Step " << stepL << "," << stepE <<std::endl;
127 
128  G4Track * theTrack = aStep->GetTrack();
129  G4int theTrackID = theTrack->GetTrackID();
130  G4double theCharge = theTrack->GetDynamicParticle()->GetCharge();
131  G4String particleType = theTrack->GetDefinition()->GetParticleName();
132  G4int pdgcode = theTrack->GetDefinition()->GetPDGEncoding();
133 
134  G4ThreeVector vert_mom = theTrack->GetVertexMomentumDirection();
135  G4double vpx = vert_mom.x();
136  G4double vpy = vert_mom.y();
137  G4double vpz = vert_mom.z();
138  double eta = 0.5 * log( (1.+vpz) / (1.-vpz) );
139  double phi = atan2(vpy,vpx);
140 
141  G4ThreeVector hitPoint = preStepPoint->GetPosition();
142  G4ThreeVector localPoint = theTrack->GetTouchable()->GetHistory()->
143  GetTopTransform().TransformPoint(hitPoint);
144 
145  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
146  int idx = touch->GetReplicaNumber(0);
147  int idLayer = -1;
148  int thePVtype = -1;
149 
150  int historyDepth = touch->GetHistoryDepth();
151 
152  if (historyDepth > 0) {
153  std::vector<int> theReplicaNumbers(historyDepth);
154  std::vector<G4VPhysicalVolume*> thePhysicalVolumes(historyDepth);
155  std::vector<G4String> thePVnames(historyDepth);
156  std::vector<G4LogicalVolume*> theLogicalVolumes(historyDepth);
157  std::vector<G4String> theLVnames(historyDepth);
158  std::vector<G4Material*> theMaterials(historyDepth);
159  std::vector<G4String> theMaterialNames(historyDepth);
160 
161  for (int jj = 0; jj < historyDepth; jj++) {
162  theReplicaNumbers[jj] = touch->GetReplicaNumber(jj);
163  thePhysicalVolumes[jj] = touch->GetVolume(jj);
164  thePVnames[jj] = thePhysicalVolumes[jj]->GetName();
165  theLogicalVolumes[jj] = thePhysicalVolumes[jj]->GetLogicalVolume();
166  theLVnames[jj] = theLogicalVolumes[jj]->GetName();
167  theMaterials[jj] = theLogicalVolumes[jj]->GetMaterial();
168  theMaterialNames[jj] = theMaterials[jj]->GetName();
169  if (verbosity >= 2)
170  std::cout << " GHD " << jj << ": " << theReplicaNumbers[jj] << ","
171  << thePVnames[jj] << "," << theLVnames[jj] << ","
172  << theMaterialNames[jj] << std::endl;
173  }
174 
175  idLayer = theReplicaNumbers[1];
176  if (thePVnames[0] == "ZDC_EMLayer")
177  thePVtype = 1;
178  else if (thePVnames[0] == "ZDC_EMAbsorber")
179  thePVtype = 2;
180  else if (thePVnames[0] == "ZDC_EMFiber")
181  thePVtype = 3;
182  else if (thePVnames[0] == "ZDC_HadLayer")
183  thePVtype = 7;
184  else if (thePVnames[0] == "ZDC_HadAbsorber")
185  thePVtype = 8;
186  else if (thePVnames[0] == "ZDC_HadFiber")
187  thePVtype = 9;
188  else if (thePVnames[0] == "ZDC_PhobosLayer")
189  thePVtype = 11;
190  else if (thePVnames[0] == "ZDC_PhobosAbsorber")
191  thePVtype = 12;
192  else if (thePVnames[0] == "ZDC_PhobosFiber")
193  thePVtype = 13;
194  else {
195  thePVtype = 0;
196  if (verbosity >= 2)
197  std::cout << " pvtype=0 hd=" << historyDepth << " " << theReplicaNumbers[0] << ","
198  << thePVnames[0] << "," << theLVnames[0] << "," << theMaterialNames[0] << std::endl;
199  }
200  }
201  else if (historyDepth == 0) {
202  int theReplicaNumber = touch->GetReplicaNumber(0);
203  G4VPhysicalVolume* thePhysicalVolume = touch->GetVolume(0);
204  G4String thePVname = thePhysicalVolume->GetName();
205  G4LogicalVolume * theLogicalVolume = thePhysicalVolume->GetLogicalVolume();
206  G4String theLVname = theLogicalVolume->GetName();
207  G4Material * theMaterial = theLogicalVolume->GetMaterial();
208  G4String theMaterialName = theMaterial->GetName();
209  if (verbosity >= 2)
210  std::cout << " hd=0 " << theReplicaNumber << ","
211  << thePVname << "," << theLVname << ","
212  << theMaterialName << std::endl;
213  }
214  else {
215  std::cout << " hd<0: hd=" << historyDepth << std::endl;
216  }
217 
218 
219  double NCherPhot = -1;
221  zdcsteparray[ntzdcs_trackid] = (float)theTrackID;
222  zdcsteparray[ntzdcs_charge] = theCharge;
223  zdcsteparray[ntzdcs_pdgcode] = (float)pdgcode;
224  zdcsteparray[ntzdcs_x] = hitPoint.x();
225  zdcsteparray[ntzdcs_y] = hitPoint.y();
226  zdcsteparray[ntzdcs_z] = hitPoint.z();
227  zdcsteparray[ntzdcs_stepl] = stepL;
228  zdcsteparray[ntzdcs_stepe] = stepE;
231  zdcsteparray[ntzdcs_vpx] = vpx;
232  zdcsteparray[ntzdcs_vpy] = vpy;
233  zdcsteparray[ntzdcs_vpz] = vpz;
234  zdcsteparray[ntzdcs_idx] = (float)idx;
235  zdcsteparray[ntzdcs_idl] = (float)idLayer;
236  zdcsteparray[ntzdcs_pvtype] = thePVtype;
237  zdcsteparray[ntzdcs_ncherphot] = NCherPhot;
239 
240  }
241 }
242 
244  //end of event
245 
246  // Look for the Hit Collection
247  std::cout << std::endl << "ZdcTest::upDate(const EndOfEvent * evt) - event #" << (*evt)()->GetEventID()
248  << std::endl << " # of aSteps followed in event = " << stepIndex << std::endl;
249 
250  // access to the G4 hit collections
251  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
252  std::cout << " accessed all HC";
253 
254  int theZDCHCid = G4SDManager::GetSDMpointer()->GetCollectionID("ZDCHITS");
255  std::cout << " - theZDCHCid = " << theZDCHCid;
256 
257  CaloG4HitCollection* theZDCHC = (CaloG4HitCollection*) allHC->GetHC(theZDCHCid);
258  std::cout << " - theZDCHC = " << theZDCHC << std::endl;
259 
260  ZdcNumberingScheme * theZdcNumScheme = new ZdcNumberingScheme(1);
261 
262  float ETot=0., SEnergy=0.;
263  int maxTime=0;
264  int fiberID=0;
265  unsigned int unsignedfiberID=0;
266  std::map<int,float,std::less<int> > energyInFibers;
267  std::map<int,float,std::less<int> > primaries;
268  float totalEnergy = 0;
269  int nentries = theZDCHC->entries();
270  std::cout << " theZDCHC has " << nentries << " entries" << std::endl;
271 
272  if (doNTzdcevent) {
273  if (nentries > 0) {
274  for (int ihit = 0; ihit < nentries; ihit++) {
275  CaloG4Hit* caloHit = (*theZDCHC)[ihit];
276  totalEnergy += caloHit->getEnergyDeposit();
277  }
278 
279  for (int ihit = 0; ihit < nentries; ihit++) {
280  CaloG4Hit* aHit = (*theZDCHC)[ihit];
281  fiberID = aHit->getUnitID();
282  unsignedfiberID = aHit->getUnitID();
283  double enEm = aHit->getEM();
284  double enHad = aHit->getHadr();
285  math::XYZPoint hitPoint = aHit->getPosition();
286  double hitEnergy = aHit->getEnergyDeposit();
287  if (verbosity >= 1)
288  std::cout << " entry #" << ihit << ": fiberID=0x" << std::hex
289  << fiberID << std::dec << "; enEm=" << enEm
290  << "; enHad=" << enHad << "; hitEnergy=" << hitEnergy
291  << "z=" << hitPoint.z() << std::endl;
292  energyInFibers[fiberID]+= enEm + enHad;
293  primaries[aHit->getTrackID()]+= enEm + enHad;
294  float time = aHit->getTimeSliceID();
295  if (time > maxTime) maxTime=(int)time;
296 
297  int thesubdet, thelayer, thefiber, thechannel, thez;
298  theZdcNumScheme->unpackZdcIndex(fiberID, thesubdet, thelayer, thefiber, thechannel, thez);
299  int unsignedsubdet, unsignedlayer, unsignedfiber, unsignedchannel, unsignedz;
300  theZdcNumScheme->unpackZdcIndex(unsignedfiberID, unsignedsubdet,
301  unsignedlayer, unsignedfiber, unsignedchannel, unsignedz);
302 
303  // unsigned int packidx1 = packZdcIndex(thesubdet, thelayer, thefiber, thechannel, thez);
304  // unsigned int packidx1 = packZdcIndex(thesubdet, thelayer, thefiber, thechannel, thez);
305  // unsigned int packidx1 = packZdcIndex(thesubdet, thelayer, thefiber, thechannel, thez);
306  // unsigned int packidx1 = packZdcIndex(thesubdet, thelayer, thefiber, thechannel, thez);
307 
309  zdceventarray[ntzdce_ihit] = (float)ihit;
310  zdceventarray[ntzdce_fiberid] = (float)fiberID;
311  zdceventarray[ntzdce_zside] = (float)thez;
312  zdceventarray[ntzdce_subdet] = (float)thesubdet;
313  zdceventarray[ntzdce_layer] = (float)thelayer;
314  zdceventarray[ntzdce_fiber] = (float)thefiber;
315  zdceventarray[ntzdce_channel] = (float)thechannel;
316  zdceventarray[ntzdce_enem] = enEm;
317  zdceventarray[ntzdce_enhad] = enHad;
318  zdceventarray[ntzdce_hitenergy] = hitEnergy;
319  zdceventarray[ntzdce_x] = hitPoint.x();
320  zdceventarray[ntzdce_y] = hitPoint.y();
321  zdceventarray[ntzdce_z] = hitPoint.z();
323  zdceventarray[ntzdce_etot] = totalEnergy;
325  }
326 
327  for (std::map<int,float,std::less<int> >::iterator is = energyInFibers.begin();
328  is!= energyInFibers.end(); is++)
329  {
330  ETot = (*is).second;
331  SEnergy += ETot;
332  }
333 
334  // Find Primary info:
335  int trackID = 0;
336  G4PrimaryParticle* thePrim=0;
337  G4int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
338  std::cout << "Event has " << nvertex << " vertex" << std::endl;
339  if (nvertex==0)
340  std::cout << "ZdcTest End Of Event ERROR: no vertex" << std::endl;
341 
342  for (int i = 0 ; i<nvertex; i++) {
343 
344  G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
345  if (avertex == 0)
346  std::cout << "ZdcTest End Of Event ERR: pointer to vertex = 0"
347  << std::endl;
348  std::cout << "Vertex number :" <<i << std::endl;
349  int npart = avertex->GetNumberOfParticle();
350  if (npart ==0)
351  std::cout << "ZdcTest End Of Event ERR: no primary!" << std::endl;
352  if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
353  }
354 
355  double px=0.,py=0.,pz=0.;
356  double pInit = 0.;
357 
358  if (thePrim != 0) {
359  px = thePrim->GetPx();
360  py = thePrim->GetPy();
361  pz = thePrim->GetPz();
362  pInit = sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
363  if (pInit==0) {
364  std::cout << "ZdcTest End Of Event ERR: primary has p=0 " << std::endl;
365  }
366  } else {
367  std::cout << "ZdcTest End Of Event ERR: could not find primary "
368  << std::endl;
369  }
370 
371  } // nentries > 0
372 
373 
374  } // createNTzdcevent
375 
376  int iEvt = (*evt)()->GetEventID();
377  if (iEvt < 10)
378  std::cout << " ZdcTest Event " << iEvt << std::endl;
379  else if ((iEvt < 100) && (iEvt%10 == 0))
380  std::cout << " ZdcTest Event " << iEvt << std::endl;
381  else if ((iEvt < 1000) && (iEvt%100 == 0))
382  std::cout << " ZdcTest Event " << iEvt << std::endl;
383  else if ((iEvt < 10000) && (iEvt%1000 == 0))
384  std::cout << " ZdcTest Event " << iEvt << std::endl;
385 }
386 
388 
390  if (doNTzdcstep) {
391  zdcOutputStepFile->cd();
392  zdcstepntuple->Write();
393  std::cout << "ZdcTestAnalysis: Ntuple step written for event: "<<eventIndex<<std::endl;
394  zdcOutputStepFile->Close();
395  std::cout << "ZdcTestAnalysis: Step file closed" << std::endl;
396  }
397 
398 
399  if (doNTzdcevent) {
400  zdcOutputEventFile->cd();
401  zdceventntuple->Write("",TObject::kOverwrite);
402  std::cout << "ZdcTestAnalysis: Ntuple event written for event: "<<eventIndex<<std::endl;
403  zdcOutputEventFile->Close();
404  std::cout << "ZdcTestAnalysis: Event file closed" << std::endl;
405  }
406 
407 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
math::XYZPoint getPosition() const
Definition: CaloG4Hit.h:56
static std::vector< std::string > checklist log
Float_t zdceventarray[16]
TFile * zdcOutputEventFile
virtual ~ZdcTestAnalysis()
double npart
Definition: HydjetWrapper.h:44
ntzdce_elements
T eta() const
std::string stepNtFileName
TFile * zdcOutputStepFile
T sqrt(T t)
Definition: SSEVec.h:48
void update(const BeginOfJob *run)
This routine will be called when the appropriate signal arrives.
TNtuple * zdceventntuple
TNtuple * zdcstepntuple
ntzdcs_elements
int getTrackID() const
Definition: CaloG4Hit.h:68
ZdcTestAnalysis(const edm::ParameterSet &p)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
static void unpackZdcIndex(const unsigned int &idx, int &subDet, int &layer, int &fiber, int &channel, int &z)
int getTimeSliceID() const
Definition: CaloG4Hit.h:71
double getEM() const
Definition: CaloG4Hit.h:59
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
std::string eventNtFileName
Float_t zdcsteparray[18]
tuple cout
Definition: gather_cfg.py:121
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
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81
Definition: DDAxes.h:10