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:
9 #include "G4SDManager.hh"
10 #include "G4Step.hh"
11 #include "G4Track.hh"
12 #include "G4VProcess.hh"
14 #include "G4ios.hh"
15 #include "G4Cerenkov.hh"
16 #include "G4ParticleTable.hh"
17 #include "CLHEP/Units/GlobalSystemOfUnits.h"
18 #include "Randomize.hh"
19 #include "G4Poisson.hh"
20 
21 ZdcSD::ZdcSD(G4String name, const DDCompactView & cpv,
23  edm::ParameterSet const & p,const SimTrackManager* manager) :
24  CaloSD(name, cpv, clg, p, manager), numberingScheme(0) {
25  edm::ParameterSet m_ZdcSD = p.getParameter<edm::ParameterSet>("ZdcSD");
26  useShowerLibrary = m_ZdcSD.getParameter<bool>("UseShowerLibrary");
27  useShowerHits = m_ZdcSD.getParameter<bool>("UseShowerHits");
28  zdcHitEnergyCut = m_ZdcSD.getParameter<double>("ZdcHitEnergyCut")*GeV;
29  verbosity = m_ZdcSD.getParameter<int>("Verbosity");
30  int verbn = verbosity/10;
31  verbosity %= 10;
32  ZdcNumberingScheme* scheme;
33  scheme = new ZdcNumberingScheme(verbn);
34  setNumberingScheme(scheme);
35 
36  edm::LogInfo("ForwardSim")
37  << "***************************************************\n"
38  << "* *\n"
39  << "* Constructing a ZdcSD with name " << name <<" *\n"
40  << "* *\n"
41  << "***************************************************";
42 
43  edm::LogInfo("ForwardSim")
44  << "\nUse of shower library is set to "
45  << useShowerLibrary
46  << "\nUse of Shower hits method is set to "
47  << useShowerHits;
48 
49  edm::LogInfo("ForwardSim")
50  << "\nEnergy Threshold Cut set to "
51  << zdcHitEnergyCut/GeV
52  <<" (GeV)";
53 
54  if(useShowerLibrary){
55  showerLibrary = new ZdcShowerLibrary(name, cpv, p);
56  }
57 }
58 
60 
62  if(showerLibrary)delete showerLibrary;
63 
64  edm::LogInfo("ForwardSim")
65  <<"end of ZdcSD\n";
66 }
67 
69  if(useShowerLibrary){
70  G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
71  showerLibrary->initRun(theParticleTable);
72  }
73  hits.clear();
74 }
75 
76 bool ZdcSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
77 
78  NaNTrap( aStep ) ;
79 
80  if (aStep == NULL) {
81  return true;
82  } else {
83  if(useShowerLibrary){
84  getFromLibrary(aStep);
85  }
86  if(useShowerHits){
87  if (getStepInfo(aStep)) {
88  if (hitExists() == false && edepositEM+edepositHAD>0.)
90  }
91  }
92  }
93  return true;
94 }
95 
96 void ZdcSD::getFromLibrary (G4Step* aStep) {
97  bool ok = true;
98 
99  preStepPoint = aStep->GetPreStepPoint();
100  theTrack = aStep->GetTrack();
101 
102  double etrack = preStepPoint->GetKineticEnergy();
103  int primaryID = setTrackID(aStep);
104 
105  hits.clear();
106 
107  /*
108  if (etrack >= zdcHitEnergyCut) {
109  primaryID = theTrack->GetTrackID();
110  } else {
111  primaryID = theTrack->GetParentID();
112  if (primaryID == 0) primaryID = theTrack->GetTrackID();
113  }
114  */
115 
116  // Reset entry point for new primary
117  posGlobal = preStepPoint->GetPosition();
118  resetForNewPrimary(posGlobal, etrack);
119 
120  if (etrack >= zdcHitEnergyCut){
121  // create hits only if above threshold
122 
123  LogDebug("ForwardSim")
124  //std::cout
125  <<"----------------New track------------------------------\n"
126  <<"Incident EnergyTrack: "<<etrack<< " MeV \n"
127  <<"Zdc Cut Energy for Hits: "<<zdcHitEnergyCut<<" MeV \n"
128  << "ZdcSD::getFromLibrary " <<hits.size() <<" hits for "
129  << GetName() << " of " << primaryID << " with "
130  << theTrack->GetDefinition()->GetParticleName() << " of "
131  << preStepPoint->GetKineticEnergy()<< " MeV\n";
132 
133  hits.swap(showerLibrary->getHits(aStep, ok));
134  }
135 
136  entrancePoint = preStepPoint->GetPosition();
137  for (unsigned int i=0; i<hits.size(); i++) {
138  posGlobal = hits[i].position;
139  entranceLocal = hits[i].entryLocal;
140  double time = hits[i].time;
141  unsigned int unitID = hits[i].detID;
142  edepositHAD = hits[i].DeHad;
143  edepositEM = hits[i].DeEM;
144  currentID.setID(unitID, time, primaryID);
145 
146  // check if it is in the same unit and timeslice as the previous on
147  if (currentID == previousID) {
149  } else {
151  }
152 
153  // currentHit->setPosition(hitPoint.x(),hitPoint.y(),hitPoint.z());
154  // currentHit->setEM(eEM);
155  // currentHit->setHadr(eHAD);
156  currentHit->setIncidentEnergy(etrack);
157  // currentHit->setEntryLocal(hitEntry.x(),hitEntry.y(),hitEntry.z());
158 
159  LogDebug("ForwardSim") << "ZdcSD: Final Hit number:"<<i<<"-->"
160  <<"New HitID: "<<currentHit->getUnitID()
161  <<" New Hit trackID: "<<currentHit->getTrackID()
162  <<" New EM Energy: "<<currentHit->getEM()/GeV
163  <<" New HAD Energy: "<<currentHit->getHadr()/GeV
164  <<" New HitEntryPoint: "<<currentHit->getEntryLocal()
165  <<" New IncidentEnergy: "<<currentHit->getIncidentEnergy()/GeV
166  <<" New HitPosition: "<<posGlobal;
167  }
168 
169  //Now kill the current track
170  if (ok) {
171  theTrack->SetTrackStatus(fStopAndKill);
172  G4TrackVector tv = *(aStep->GetSecondary());
173  for (unsigned int kk=0; kk<tv.size(); kk++) {
174  if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
175  tv[kk]->SetTrackStatus(fStopAndKill);
176  }
177  }
178 }
179 
180 double ZdcSD::getEnergyDeposit(G4Step * aStep, edm::ParameterSet const & p ) {
181 
182  float NCherPhot = 0.;
183  //std::cout<<"I go through here"<<std::endl;
184 
185  if (aStep == NULL) {
186  LogDebug("ForwardSim") << "ZdcSD:: getEnergyDeposit: aStep is NULL!";
187  return 0;
188  } else {
189  // preStepPoint information
190  G4SteppingControl stepControlFlag = aStep->GetControlFlag();
191  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
192  G4VPhysicalVolume* currentPV = preStepPoint->GetPhysicalVolume();
193  G4String nameVolume = currentPV->GetName();
194 
195  G4ThreeVector hitPoint = preStepPoint->GetPosition();
196  G4ThreeVector hit_mom = preStepPoint->GetMomentumDirection();
197  G4double stepL = aStep->GetStepLength()/cm;
198  G4double beta = preStepPoint->GetBeta();
199  G4double charge = preStepPoint->GetCharge();
200 
201  // G4VProcess* curprocess = preStepPoint->GetProcessDefinedStep();
202  // G4String namePr = preStepPoint->GetProcessDefinedStep()->GetProcessName();
203  // G4LogicalVolume* lv = currentPV->GetLogicalVolume();
204  // G4Material* mat = lv->GetMaterial();
205  // G4double rad = mat->GetRadlen();
206 
207  // postStepPoint information
208  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
209  G4VPhysicalVolume* postPV = postStepPoint->GetPhysicalVolume();
210  G4String postnameVolume = postPV->GetName();
211 
212  // theTrack information
213  G4Track* theTrack = aStep->GetTrack();
214  G4String particleType = theTrack->GetDefinition()->GetParticleName();
215  G4int primaryID = theTrack->GetTrackID();
216  G4double entot = theTrack->GetTotalEnergy();
217  G4ThreeVector vert_mom = theTrack->GetVertexMomentumDirection();
218  G4ThreeVector localPoint = theTrack->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
219 
220  // calculations
221  float costheta = vert_mom.z()/sqrt(vert_mom.x()*vert_mom.x()+
222  vert_mom.y()*vert_mom.y()+
223  vert_mom.z()*vert_mom.z());
224  float theta = acos(std::min(std::max(costheta,float(-1.)),float(1.)));
225  float eta = -log(tan(theta/2));
226  float phi = -100.;
227  if (vert_mom.x() != 0) phi = atan2(vert_mom.y(),vert_mom.x());
228  if (phi < 0.) phi += twopi;
229 
230  // Get the total energy deposit
231  double stepE = aStep->GetTotalEnergyDeposit();
232  LogDebug("ForwardSim")
233  << "ZdcSD:: getEnergyDeposit: "
234  <<"*****************HHHHHHHHHHHHHHHHHHHHHHHHHHLLLLLLLLLlllllllllll&&&&&&&&&&\n"
235  << " preStepPoint: " << nameVolume << "," << stepL << "," << stepE
236  << "," << beta << "," << charge << "\n"
237  << " postStepPoint: " << postnameVolume << "," << costheta << ","
238  << theta << "," << eta << "," << phi << "," << particleType << ","
239  << primaryID;
240 
241  float bThreshold = 0.67;
242  int status = 0;
243  if ((beta > bThreshold) && (charge != 0) && (nameVolume == "ZDC_EMFiber" || nameVolume == "ZDC_HadFiber")) {
244  status = 1;
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  status = 0;
397  if (beta <= bThreshold)
398  LogDebug("ForwardSim")
399  << "ZdcSD:: getEnergyDeposit: fail beta=" << beta;
400  if (charge == 0)
401  LogDebug("ForwardSim")
402  << "ZdcSD:: getEnergyDeposit: fail charge=0";
403  if ( !(nameVolume == "ZDC_EMFiber" || nameVolume == "ZDC_HadFiber") )
404  LogDebug("ForwardSim")
405  << "ZdcSD:: getEnergyDeposit: fail nv=" << nameVolume;
406  }
407 
408  return NCherPhot;
409  }
410 }
411 
412 uint32_t ZdcSD::setDetUnitId(G4Step* aStep) {
413  uint32_t returnNumber = 0;
414  if(numberingScheme != 0)returnNumber = numberingScheme->getUnitID(aStep);
415  // edm: return (numberingScheme == 0 ? 0 : numberingScheme->getUnitID(aStep));
416  return returnNumber;
417 }
418 
420  if (scheme != 0) {
421  edm::LogInfo("ForwardSim") << "ZdcSD: updates numbering scheme for "
422  << GetName();
423  if (numberingScheme) delete numberingScheme;
424  numberingScheme = scheme;
425  }
426 }
427 
428 int ZdcSD::setTrackID (G4Step* aStep) {
429  theTrack = aStep->GetTrack();
430  double etrack = preStepPoint->GetKineticEnergy();
431  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
432  int primaryID = trkInfo->getIDonCaloSurface();
433  if (primaryID == 0) {
434 #ifdef DebugLog
435  LogDebug("ZdcSD") << "ZdcSD: Problem with primaryID **** set by force "
436  << "to TkID **** " << theTrack->GetTrackID();
437 #endif
438  primaryID = theTrack->GetTrackID();
439  }
440  if (primaryID != previousID.trackID())
441  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
442  return primaryID;
443 }
#define LogDebug(id)
float edepositEM
Definition: CaloSD.h:119
const double beta
void setNumberingScheme(ZdcNumberingScheme *scheme)
Definition: ZdcSD.cc:419
double thFibDir
Definition: ZdcSD.h:36
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void updateHit(CaloG4Hit *)
Definition: CaloSD.cc:453
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:59
virtual double getEnergyDeposit(G4Step *, edm::ParameterSet const &)
Definition: ZdcSD.cc:180
ZdcNumberingScheme * numberingScheme
Definition: ZdcSD.h:39
double getIncidentEnergy() const
Definition: CaloG4Hit.h:65
#define NULL
Definition: scimark2.h:8
Geom::Theta< T > theta() const
#define min(a, b)
Definition: mlp_lapack.h:161
G4ThreeVector posGlobal
Definition: CaloSD.h:111
bool useShowerHits
Definition: ZdcSD.h:34
std::vector< Hit > & getHits(G4Step *aStep, bool &ok)
type of data representation of DDCompactView
Definition: DDCompactView.h:81
double charge(const std::vector< uint8_t > &Ampls)
T eta() const
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
const T & max(const T &a, const T &b)
float edepositHAD
Definition: CaloSD.h:119
T sqrt(T t)
Definition: SSEVec.h:28
int trackID() const
Definition: CaloHitID.h:25
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
CaloHitID previousID
Definition: CaloSD.h:115
CaloG4Hit * currentHit
Definition: CaloSD.h:126
void NaNTrap(G4Step *step)
virtual void initRun()
Definition: ZdcSD.cc:68
G4Track * theTrack
Definition: CaloSD.h:116
virtual G4bool getStepInfo(G4Step *aStep)
Definition: CaloSD.cc:267
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:44
G4StepPoint * preStepPoint
Definition: CaloSD.h:118
void getFromLibrary(G4Step *step)
Definition: ZdcSD.cc:96
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:76
Log< T >::type log(const T &t)
Definition: Log.h:22
CaloHitID currentID
Definition: CaloSD.h:115
double getEM() const
Definition: CaloG4Hit.h:59
double a
Definition: hdecay.h:121
double zdcHitEnergyCut
Definition: ZdcSD.h:37
G4bool hitExists()
Definition: CaloSD.cc:319
math::XYZPoint getEntryLocal() const
Definition: CaloG4Hit.h:53
int setTrackID(G4Step *step)
Definition: ZdcSD.cc:428
void resetForNewPrimary(G4ThreeVector, double)
Definition: CaloSD.cc:468
ZdcSD(G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: ZdcSD.cc:21
double pi
uint32_t getUnitID() const
Definition: CaloG4Hit.h:69
tuple status
Definition: ntuplemaker.py:245
void initRun(G4ParticleTable *theParticleTable)
virtual uint32_t setDetUnitId(G4Step *step)
Definition: ZdcSD.cc:412
bool useShowerLibrary
Definition: ZdcSD.h:34
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:375
G4ThreeVector entrancePoint
Definition: CaloSD.h:109
G4ThreeVector entranceLocal
Definition: CaloSD.h:110
double getHadr() const
Definition: CaloG4Hit.h:62
ZdcShowerLibrary * showerLibrary
Definition: ZdcSD.h:38
Definition: DDAxes.h:10