CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZdcSD.cc
Go to the documentation of this file.
1 // File: ZdcSD.cc
3 // Date: 03.01
4 // Description: Sensitive Detector class for Zdc
5 // Modifications:
10 
11 #include "G4SDManager.hh"
12 #include "G4Step.hh"
13 #include "G4Track.hh"
14 #include "G4VProcess.hh"
15 #include "G4ios.hh"
16 #include "G4Cerenkov.hh"
17 #include "G4ParticleTable.hh"
18 #include "CLHEP/Units/GlobalSystemOfUnits.h"
19 #include "CLHEP/Units/GlobalPhysicalConstants.h"
20 #include "Randomize.hh"
21 #include "G4Poisson.hh"
22 
23 ZdcSD::ZdcSD(G4String name, const DDCompactView & cpv,
24  const SensitiveDetectorCatalog & clg,
25  edm::ParameterSet const & p,const SimTrackManager* manager) :
26  CaloSD(name, cpv, clg, p, manager), numberingScheme(0) {
27  edm::ParameterSet m_ZdcSD = p.getParameter<edm::ParameterSet>("ZdcSD");
28  useShowerLibrary = m_ZdcSD.getParameter<bool>("UseShowerLibrary");
29  useShowerHits = m_ZdcSD.getParameter<bool>("UseShowerHits");
30  zdcHitEnergyCut = m_ZdcSD.getParameter<double>("ZdcHitEnergyCut")*GeV;
31  verbosity = m_ZdcSD.getParameter<int>("Verbosity");
32  int verbn = verbosity/10;
33  verbosity %= 10;
34  ZdcNumberingScheme* scheme;
35  scheme = new ZdcNumberingScheme(verbn);
36  setNumberingScheme(scheme);
37 
38  edm::LogInfo("ForwardSim")
39  << "***************************************************\n"
40  << "* *\n"
41  << "* Constructing a ZdcSD with name " << name <<" *\n"
42  << "* *\n"
43  << "***************************************************";
44 
45  edm::LogInfo("ForwardSim")
46  << "\nUse of shower library is set to "
47  << useShowerLibrary
48  << "\nUse of Shower hits method is set to "
49  << useShowerHits;
50 
51  edm::LogInfo("ForwardSim")
52  << "\nEnergy Threshold Cut set to "
53  << zdcHitEnergyCut/GeV
54  <<" (GeV)";
55 
56  if(useShowerLibrary){
57  showerLibrary = new ZdcShowerLibrary(name, cpv, p);
58  }
59 }
60 
62 
64  if(showerLibrary)delete showerLibrary;
65 
66  edm::LogInfo("ForwardSim")
67  <<"end of ZdcSD\n";
68 }
69 
71  if(useShowerLibrary){
72  G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
73  showerLibrary->initRun(theParticleTable);
74  }
75  hits.clear();
76 }
77 
78 bool ZdcSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
79 
80  NaNTrap( aStep ) ;
81 
82  if (aStep == NULL) {
83  return true;
84  } else {
85  if(useShowerLibrary){
86  getFromLibrary(aStep);
87  }
88  if(useShowerHits){
89  if (getStepInfo(aStep)) {
90  if (hitExists() == false && edepositEM+edepositHAD>0.)
92  }
93  }
94  }
95  return true;
96 }
97 
98 void ZdcSD::getFromLibrary (G4Step* aStep) {
99  bool ok = true;
100 
101  preStepPoint = aStep->GetPreStepPoint();
102  theTrack = aStep->GetTrack();
103 
104  double etrack = preStepPoint->GetKineticEnergy();
105  int primaryID = setTrackID(aStep);
106 
107  hits.clear();
108 
109  /*
110  if (etrack >= zdcHitEnergyCut) {
111  primaryID = theTrack->GetTrackID();
112  } else {
113  primaryID = theTrack->GetParentID();
114  if (primaryID == 0) primaryID = theTrack->GetTrackID();
115  }
116  */
117 
118  // Reset entry point for new primary
119  posGlobal = preStepPoint->GetPosition();
120  resetForNewPrimary(posGlobal, etrack);
121 
122  if (etrack >= zdcHitEnergyCut){
123  // create hits only if above threshold
124 
125  LogDebug("ForwardSim")
126  //std::cout
127  <<"----------------New track------------------------------\n"
128  <<"Incident EnergyTrack: "<<etrack<< " MeV \n"
129  <<"Zdc Cut Energy for Hits: "<<zdcHitEnergyCut<<" MeV \n"
130  << "ZdcSD::getFromLibrary " <<hits.size() <<" hits for "
131  << GetName() << " of " << primaryID << " with "
132  << theTrack->GetDefinition()->GetParticleName() << " of "
133  << preStepPoint->GetKineticEnergy()<< " MeV\n";
134 
135  hits.swap(showerLibrary->getHits(aStep, ok));
136  }
137 
138  entrancePoint = preStepPoint->GetPosition();
139  for (unsigned int i=0; i<hits.size(); i++) {
140  posGlobal = hits[i].position;
141  entranceLocal = hits[i].entryLocal;
142  double time = hits[i].time;
143  unsigned int unitID = hits[i].detID;
144  edepositHAD = hits[i].DeHad;
145  edepositEM = hits[i].DeEM;
146  currentID.setID(unitID, time, primaryID);
147 
148  // check if it is in the same unit and timeslice as the previous on
149  if (currentID == previousID) {
151  } else {
153  }
154 
155  // currentHit->setPosition(hitPoint.x(),hitPoint.y(),hitPoint.z());
156  // currentHit->setEM(eEM);
157  // currentHit->setHadr(eHAD);
158  currentHit->setIncidentEnergy(etrack);
159  // currentHit->setEntryLocal(hitEntry.x(),hitEntry.y(),hitEntry.z());
160 
161  LogDebug("ForwardSim") << "ZdcSD: Final Hit number:"<<i<<"-->"
162  <<"New HitID: "<<currentHit->getUnitID()
163  <<" New Hit trackID: "<<currentHit->getTrackID()
164  <<" New EM Energy: "<<currentHit->getEM()/GeV
165  <<" New HAD Energy: "<<currentHit->getHadr()/GeV
166  <<" New HitEntryPoint: "<<currentHit->getEntryLocal()
167  <<" New IncidentEnergy: "<<currentHit->getIncidentEnergy()/GeV
168  <<" New HitPosition: "<<posGlobal;
169  }
170 
171  //Now kill the current track
172  if (ok) {
173  theTrack->SetTrackStatus(fStopAndKill);
174  G4TrackVector tv = *(aStep->GetSecondary());
175  for (unsigned int kk=0; kk<tv.size(); kk++) {
176  if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
177  tv[kk]->SetTrackStatus(fStopAndKill);
178  }
179  }
180 }
181 
182 double ZdcSD::getEnergyDeposit(G4Step * aStep, edm::ParameterSet const & p ) {
183 
184  float NCherPhot = 0.;
185  //std::cout<<"I go through here"<<std::endl;
186 
187  if (aStep == NULL) {
188  LogDebug("ForwardSim") << "ZdcSD:: getEnergyDeposit: aStep is NULL!";
189  return 0;
190  } else {
191  // preStepPoint information
192  G4SteppingControl stepControlFlag = aStep->GetControlFlag();
193  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
194  G4VPhysicalVolume* currentPV = preStepPoint->GetPhysicalVolume();
195  G4String nameVolume = currentPV->GetName();
196 
197  G4ThreeVector hitPoint = preStepPoint->GetPosition();
198  G4ThreeVector hit_mom = preStepPoint->GetMomentumDirection();
199  G4double stepL = aStep->GetStepLength()/cm;
200  G4double beta = preStepPoint->GetBeta();
201  G4double charge = preStepPoint->GetCharge();
202 
203  // G4VProcess* curprocess = preStepPoint->GetProcessDefinedStep();
204  // G4String namePr = preStepPoint->GetProcessDefinedStep()->GetProcessName();
205  // G4LogicalVolume* lv = currentPV->GetLogicalVolume();
206  // G4Material* mat = lv->GetMaterial();
207  // G4double rad = mat->GetRadlen();
208 
209  // postStepPoint information
210  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
211  G4VPhysicalVolume* postPV = postStepPoint->GetPhysicalVolume();
212  G4String postnameVolume = postPV->GetName();
213 
214  // theTrack information
215  G4Track* theTrack = aStep->GetTrack();
216  G4String particleType = theTrack->GetDefinition()->GetParticleName();
217  G4int primaryID = theTrack->GetTrackID();
218  G4double entot = theTrack->GetTotalEnergy();
219  G4ThreeVector vert_mom = theTrack->GetVertexMomentumDirection();
220  G4ThreeVector localPoint = theTrack->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
221 
222  // calculations
223  float costheta = vert_mom.z()/sqrt(vert_mom.x()*vert_mom.x()+
224  vert_mom.y()*vert_mom.y()+
225  vert_mom.z()*vert_mom.z());
226  float theta = acos(std::min(std::max(costheta,float(-1.)),float(1.)));
227  float eta = -log(tan(theta/2));
228  float phi = -100.;
229  if (vert_mom.x() != 0) phi = atan2(vert_mom.y(),vert_mom.x());
230  if (phi < 0.) phi += twopi;
231 
232  // Get the total energy deposit
233  double stepE = aStep->GetTotalEnergyDeposit();
234  LogDebug("ForwardSim")
235  << "ZdcSD:: getEnergyDeposit: "
236  <<"*****************HHHHHHHHHHHHHHHHHHHHHHHHHHLLLLLLLLLlllllllllll&&&&&&&&&&\n"
237  << " preStepPoint: " << nameVolume << "," << stepL << "," << stepE
238  << "," << beta << "," << charge << "\n"
239  << " postStepPoint: " << postnameVolume << "," << costheta << ","
240  << theta << "," << eta << "," << phi << "," << particleType << ","
241  << primaryID;
242 
243  float bThreshold = 0.67;
244  if ((beta > bThreshold) && (charge != 0) && (nameVolume == "ZDC_EMFiber" || nameVolume == "ZDC_HadFiber")) {
245  LogDebug("ForwardSim") << "ZdcSD:: getEnergyDeposit: pass ";
246 
247  float nMedium = 1.4925;
248  // float photEnSpectrDL = 10714.285714;
249  // photEnSpectrDL = (1./400.nm-1./700.nm)*10000000.cm/nm; /* cm-1 */
250 
251  float photEnSpectrDE = 1.24;
252  // E = 2pi*(1./137.)*(eV*cm/370.)/lambda = 12.389184*(eV*cm)/lambda
253  // Emax = 12.389184*(eV*cm)/400nm*10-7cm/nm = 3.01 eV
254  // Emin = 12.389184*(eV*cm)/700nm*10-7cm/nm = 1.77 eV
255  // delE = Emax - Emin = 1.24 eV
256 
257  float effPMTandTransport = 0.15;
258 
259  // Check these values
260  float thFullRefl = 23.;
261  float thFullReflRad = thFullRefl*pi/180.;
262 
263  edm::ParameterSet m_ZdcSD = p.getParameter<edm::ParameterSet>("ZdcSD");
264  thFibDir = m_ZdcSD.getParameter<double>("FiberDirection");
265  //float thFibDir = 90.;
266  float thFibDirRad = thFibDir*pi/180.;
267 
268  // at which theta the point is located:
269  // float th1 = hitPoint.theta();
270 
271  // theta of charged particle in LabRF(hit momentum direction):
272  float costh = hit_mom.z()/sqrt(hit_mom.x()*hit_mom.x()+
273  hit_mom.y()*hit_mom.y()+
274  hit_mom.z()*hit_mom.z());
275  float th = acos(std::min(std::max(costh,float(-1.)),float(1.)));
276  // just in case (can do both standard ranges of phi):
277  if (th < 0.) th += twopi;
278 
279  // theta of cone with Cherenkov photons w.r.t.direction of charged part.:
280  float costhcher =1./(nMedium*beta);
281  float thcher = acos(std::min(std::max(costhcher,float(-1.)),float(1.)));
282 
283  // diff thetas of charged part. and quartz direction in LabRF:
284  float DelFibPart = fabs(th - thFibDirRad);
285 
286  // define real distances:
287  float d = fabs(tan(th)-tan(thFibDirRad));
288 
289  // float a = fabs(tan(thFibDirRad)-tan(thFibDirRad+thFullReflRad));
290  // float r = fabs(tan(th)-tan(th+thcher));
291  float a = tan(thFibDirRad)+tan(fabs(thFibDirRad-thFullReflRad));
292  float r = tan(th)+tan(fabs(th-thcher));
293 
294  // std::cout.testOut << " d=|tan(" << th << ")-tan(" << thFibDirRad << ")| "
295  // << "=|" << tan(th) << "-" << tan(thFibDirRad) << "| = " << d;
296  // std::cout.testOut << " a=tan(" << thFibDirRad << ")=" << tan(thFibDirRad)
297  // << " + tan(|" << thFibDirRad << " - " << thFullReflRad << "|)="
298  // << tan(fabs(thFibDirRad-thFullReflRad)) << " = " << a;
299  // std::cout.testOut << " r=tan(" << th << ")=" << tan(th) << " + tan(|" << th
300  // << " - " << thcher << "|)=" << tan(fabs(th-thcher)) << " = " << r;
301 
302  // define losses d_qz in cone of full reflection inside quartz direction
303  float d_qz = -1;
304  float variant = -1;
305 
306  // if (d > (r+a))
307  if (DelFibPart > (thFullReflRad + thcher) ) {
308  variant = 0.; d_qz = 0.;
309  } else {
310  // if ((DelFibPart + thcher) < thFullReflRad ) [(d+r) < a]
311  if ((th + thcher) < (thFibDirRad+thFullReflRad) && (th - thcher) > (thFibDirRad-thFullReflRad) ) {
312  variant = 1.; d_qz = 1.;
313  } else {
314  // if ((thcher - DelFibPart ) > thFullReflRad ) [(r-d) > a]
315  if ((thFibDirRad + thFullReflRad) < (th + thcher) && (thFibDirRad - thFullReflRad) > (th - thcher) ) {
316  variant = 2.; d_qz = 0.;
317  } else {
318  // if ((thcher + DelFibPart ) > thFullReflRad && thcher < (DelFibPart+thFullReflRad) ) { [(r+d) > a && (r-d) < a)]
319  variant = 3.; // d_qz is calculated below
320 
321  // use crossed length of circles(cone projection) - dC1/dC2 :
322  float arg_arcos = 0.;
323  float tan_arcos = 2.*a*d;
324  if (tan_arcos != 0.) arg_arcos =(r*r-a*a-d*d)/tan_arcos;
325  // std::cout.testOut << " d_qz: " << r << "," << a << "," << d << " " << tan_arcos << " " << arg_arcos;
326  arg_arcos = fabs(arg_arcos);
327  // std::cout.testOut << "," << arg_arcos;
328  float th_arcos = acos(std::min(std::max(arg_arcos,float(-1.)),float(1.)));
329  // std::cout.testOut << " " << th_arcos;
330  d_qz = th_arcos/pi/2.;
331  // std::cout.testOut << " " << d_qz;
332  d_qz = fabs(d_qz);
333  // std::cout.testOut << "," << d_qz;
334  }
335  }
336  }
337 
338  // std::cout<< std::endl;
339 
340  double meanNCherPhot = 0.;
341  G4int poissNCherPhot = 0;
342  if (d_qz > 0) {
343  meanNCherPhot = 370.*charge*charge*( 1. - 1./(nMedium*nMedium*beta*beta) ) * photEnSpectrDE * stepL;
344 
345  // dLamdX: meanNCherPhot = (2.*pi/137.)*charge*charge*
346  // ( 1. - 1./(nMedium*nMedium*beta*beta) ) * photEnSpectrDL * stepL;
347  poissNCherPhot = (G4int) G4Poisson(meanNCherPhot);
348 
349  if (poissNCherPhot < 0) poissNCherPhot = 0;
350 
351  // NCherPhot = meanNCherPhot;
352  NCherPhot = poissNCherPhot * effPMTandTransport * d_qz;
353  }
354 
355  LogDebug("ForwardSim")
356  << "ZdcSD:: getEnergyDeposit: gED: "
357  << stepE
358  << "," << costh
359  << "," << th
360  << "," << costhcher
361  << "," << thcher
362  << "," << DelFibPart
363  << "," << d
364  << "," << a
365  << "," << r
366  << "," << hitPoint
367  << "," << hit_mom
368  << "," << stepControlFlag
369  << "," << entot
370  << "," << vert_mom
371  << "," << localPoint
372  << "," << charge
373  << "," << beta
374  << "," << stepL
375  << "," << d_qz
376  << "," << variant
377  << "," << meanNCherPhot
378  << "," << poissNCherPhot
379  << "," << NCherPhot;
380  // --constants-----------------
381  // << "," << photEnSpectrDE
382  // << "," << nMedium
383  // << "," << bThreshold
384  // << "," << thFibDirRad
385  // << "," << thFullReflRad
386  // << "," << effPMTandTransport
387  // --other variables-----------
388  // << "," << curprocess
389  // << "," << nameProcess
390  // << "," << name
391  // << "," << rad
392  // << "," << mat
393 
394  } else {
395  // determine failure mode: beta, charge, and/or nameVolume
396  if (beta <= bThreshold)
397  LogDebug("ForwardSim")
398  << "ZdcSD:: getEnergyDeposit: fail beta=" << beta;
399  if (charge == 0)
400  LogDebug("ForwardSim")
401  << "ZdcSD:: getEnergyDeposit: fail charge=0";
402  if ( !(nameVolume == "ZDC_EMFiber" || nameVolume == "ZDC_HadFiber") )
403  LogDebug("ForwardSim")
404  << "ZdcSD:: getEnergyDeposit: fail nv=" << nameVolume;
405  }
406 
407  return NCherPhot;
408  }
409 }
410 
411 uint32_t ZdcSD::setDetUnitId(G4Step* aStep) {
412  uint32_t returnNumber = 0;
413  if(numberingScheme != 0)returnNumber = numberingScheme->getUnitID(aStep);
414  // edm: return (numberingScheme == 0 ? 0 : numberingScheme->getUnitID(aStep));
415  return returnNumber;
416 }
417 
419  if (scheme != 0) {
420  edm::LogInfo("ForwardSim") << "ZdcSD: updates numbering scheme for "
421  << GetName();
422  if (numberingScheme) delete numberingScheme;
423  numberingScheme = scheme;
424  }
425 }
426 
427 int ZdcSD::setTrackID (G4Step* aStep) {
428  theTrack = aStep->GetTrack();
429  double etrack = preStepPoint->GetKineticEnergy();
430  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
431  int primaryID = trkInfo->getIDonCaloSurface();
432  if (primaryID == 0) {
433 #ifdef DebugLog
434  LogDebug("ZdcSD") << "ZdcSD: Problem with primaryID **** set by force "
435  << "to TkID **** " << theTrack->GetTrackID();
436 #endif
437  primaryID = theTrack->GetTrackID();
438  }
439  if (primaryID != previousID.trackID())
440  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
441  return primaryID;
442 }
#define LogDebug(id)
float edepositEM
Definition: CaloSD.h:121
const double beta
void setNumberingScheme(ZdcNumberingScheme *scheme)
Definition: ZdcSD.cc:418
double thFibDir
Definition: ZdcSD.h:36
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const double GeV
Definition: MathUtil.h:16
void updateHit(CaloG4Hit *)
Definition: CaloSD.cc:440
int getIDonCaloSurface() const
Definition: CaloSD.h:42
virtual unsigned int getUnitID(const G4Step *aStep) const
int verbosity
Definition: ZdcSD.h:33
std::vector< ZdcShowerLibrary::Hit > hits
Definition: ZdcSD.h:41
void setIncidentEnergy(double e)
Definition: CaloG4Hit.h:66
virtual ~ZdcSD()
Definition: ZdcSD.cc:61
virtual double getEnergyDeposit(G4Step *, edm::ParameterSet const &)
Definition: ZdcSD.cc:182
Geom::Theta< T > theta() const
ZdcNumberingScheme * numberingScheme
Definition: ZdcSD.h:39
double getIncidentEnergy() const
Definition: CaloG4Hit.h:65
#define NULL
Definition: scimark2.h:8
G4ThreeVector posGlobal
Definition: CaloSD.h:113
bool useShowerHits
Definition: ZdcSD.h:34
ZdcSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: ZdcSD.cc:23
T eta() const
std::vector< Hit > & getHits(G4Step *aStep, bool &ok)
type of data representation of DDCompactView
Definition: DDCompactView.h:77
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
const Double_t pi
tuple d
Definition: ztail.py:151
void resetForNewPrimary(const G4ThreeVector &, double)
Definition: CaloSD.cc:454
float edepositHAD
Definition: CaloSD.h:121
T sqrt(T t)
Definition: SSEVec.h:48
int trackID() const
Definition: CaloHitID.h:25
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
T min(T a, T b)
Definition: MathUtil.h:58
CaloHitID previousID
Definition: CaloSD.h:117
CaloG4Hit * currentHit
Definition: CaloSD.h:128
void NaNTrap(G4Step *step)
virtual void initRun()
Definition: ZdcSD.cc:70
G4Track * theTrack
Definition: CaloSD.h:118
virtual G4bool getStepInfo(G4Step *aStep)
Definition: CaloSD.cc:250
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:44
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
void getFromLibrary(G4Step *step)
Definition: ZdcSD.cc:98
static const G4LogicalVolume * GetVolume(const std::string &name)
int getTrackID() const
Definition: CaloG4Hit.h:68
virtual bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory)
Definition: ZdcSD.cc:78
CaloHitID currentID
Definition: CaloSD.h:117
double getEM() const
Definition: CaloG4Hit.h:59
double a
Definition: hdecay.h:121
double zdcHitEnergyCut
Definition: ZdcSD.h:37
G4bool hitExists()
Definition: CaloSD.cc:310
math::XYZPoint getEntryLocal() const
Definition: CaloG4Hit.h:53
int setTrackID(G4Step *step)
Definition: ZdcSD.cc:427
uint32_t getUnitID() const
Definition: CaloG4Hit.h:69
void initRun(G4ParticleTable *theParticleTable)
virtual uint32_t setDetUnitId(G4Step *step)
Definition: ZdcSD.cc:411
bool useShowerLibrary
Definition: ZdcSD.h:34
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:363
G4ThreeVector entrancePoint
Definition: CaloSD.h:111
G4ThreeVector entranceLocal
Definition: CaloSD.h:112
double getHadr() const
Definition: CaloG4Hit.h:62
tuple log
Definition: cmsBatch.py:347
ZdcShowerLibrary * showerLibrary
Definition: ZdcSD.h:38
Definition: DDAxes.h:10