CMS 3D CMS Logo

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:
7 #include <memory>
8 
15 
16 #include "G4SDManager.hh"
17 #include "G4Step.hh"
18 #include "G4Track.hh"
19 #include "G4VProcess.hh"
20 #include "G4ios.hh"
21 #include "G4Cerenkov.hh"
22 #include "G4ParticleTable.hh"
23 #include "CLHEP/Units/GlobalSystemOfUnits.h"
24 #include "CLHEP/Units/GlobalPhysicalConstants.h"
25 #include "Randomize.hh"
26 #include "G4Poisson.hh"
27 
28 //#define EDM_ML_DEBUG
29 
31  const SensitiveDetectorCatalog& clg,
32  edm::ParameterSet const& p,
33  const SimTrackManager* manager)
34  : CaloSD(name, clg, p, manager) {
35  edm::ParameterSet m_ZdcSD = p.getParameter<edm::ParameterSet>("ZdcSD");
36  useShowerLibrary = m_ZdcSD.getParameter<bool>("UseShowerLibrary");
37  useShowerHits = m_ZdcSD.getParameter<bool>("UseShowerHits");
38  zdcHitEnergyCut = m_ZdcSD.getParameter<double>("ZdcHitEnergyCut") * GeV;
39  thFibDir = m_ZdcSD.getParameter<double>("FiberDirection");
40  verbosity = m_ZdcSD.getParameter<int>("Verbosity");
41  int verbn = verbosity / 10;
42  verbosity %= 10;
44 
45  edm::LogVerbatim("ForwardSim") << "***************************************************\n"
46  << "* *\n"
47  << "* Constructing a ZdcSD with name " << name << " *\n"
48  << "* *\n"
49  << "***************************************************";
50 
51  edm::LogVerbatim("ForwardSim") << "\nUse of shower library is set to " << useShowerLibrary
52  << "\nUse of Shower hits method is set to " << useShowerHits;
53 
54  edm::LogVerbatim("ForwardSim") << "\nEnergy Threshold Cut set to " << zdcHitEnergyCut / GeV << " (GeV)";
55 
56  if (useShowerLibrary) {
57  showerLibrary = std::make_unique<ZdcShowerLibrary>(name, p);
58  setParameterized(true);
59  } else {
60  showerLibrary.reset(nullptr);
61  }
62 }
63 
64 void ZdcSD::initRun() { hits.clear(); }
65 
66 bool ZdcSD::getFromLibrary(const G4Step* aStep) {
67  bool ok = true;
68 
69  auto const preStepPoint = aStep->GetPreStepPoint();
70 
71  double etrack = preStepPoint->GetKineticEnergy();
72  int primaryID = setTrackID(aStep);
73 
74  hits.clear();
75 
76  // Reset entry point for new primary
77  resetForNewPrimary(aStep);
78 
79  if (etrack >= zdcHitEnergyCut) {
80  // create hits only if above threshold
81 
82 #ifdef EDM_ML_DEBUG
83  auto const theTrack = aStep->GetTrack();
84  edm::LogVerbatim("ForwardSim") << "----------------New track------------------------------\n"
85  << "Incident EnergyTrack: " << etrack << " MeV \n"
86  << "Zdc Cut Energy for Hits: " << zdcHitEnergyCut << " MeV \n"
87  << "ZdcSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
88  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
89  << etrack << " MeV\n";
90 #endif
91  hits.swap(showerLibrary.get()->getHits(aStep, ok));
92  }
93 
94  incidentEnergy = etrack;
95  entrancePoint = preStepPoint->GetPosition();
96  for (unsigned int i = 0; i < hits.size(); i++) {
97  posGlobal = hits[i].position;
98  entranceLocal = hits[i].entryLocal;
99  double time = hits[i].time;
100  unsigned int unitID = hits[i].detID;
101  edepositHAD = hits[i].DeHad;
102  edepositEM = hits[i].DeEM;
103  currentID.setID(unitID, time, primaryID, 0);
104  processHit(aStep);
105 
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("ForwardSim") << "ZdcSD: Final Hit number:" << i << "-->"
108  << "New HitID: " << currentHit->getUnitID()
109  << " New Hit trackID: " << currentHit->getTrackID()
110  << " New EM Energy: " << currentHit->getEM() / GeV
111  << " New HAD Energy: " << currentHit->getHadr() / GeV
112  << " New HitEntryPoint: " << currentHit->getEntryLocal()
113  << " New IncidentEnergy: " << currentHit->getIncidentEnergy() / GeV
114  << " New HitPosition: " << posGlobal;
115 #endif
116  }
117  return ok;
118 }
119 
120 double ZdcSD::getEnergyDeposit(const G4Step* aStep) {
121  double NCherPhot = 0.;
122 
123  // preStepPoint information
124  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
125 
126  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
127  const G4ThreeVector& hit_mom = preStepPoint->GetMomentumDirection();
128  G4double stepL = aStep->GetStepLength() / cm;
129  G4double beta = preStepPoint->GetBeta();
130  G4double charge = preStepPoint->GetCharge();
131  if (charge == 0.0)
132  return 0.0;
133 
134  // theTrack information
135  G4Track* theTrack = aStep->GetTrack();
136  G4ThreeVector localPoint = theTrack->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
137 
138 #ifdef EDM_ML_DEBUG
139  const G4ThreeVector& vert_mom = theTrack->GetVertexMomentumDirection();
140 
141  // calculations
142  float costheta =
143  vert_mom.z() / sqrt(vert_mom.x() * vert_mom.x() + vert_mom.y() * vert_mom.y() + vert_mom.z() * vert_mom.z());
144  float theta = std::acos(std::min(std::max(costheta, -1.f), 1.f));
145  float eta = -std::log(std::tan(theta * 0.5f));
146  float phi = -100.;
147  if (vert_mom.x() != 0)
148  phi = std::atan2(vert_mom.y(), vert_mom.x());
149  if (phi < 0.)
150  phi += twopi;
151 
152  // Get the total energy deposit
153  double stepE = aStep->GetTotalEnergyDeposit();
154 
155  // postStepPoint information
156  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
157  G4VPhysicalVolume* postPV = postStepPoint->GetPhysicalVolume();
158  G4VPhysicalVolume* currentPV = preStepPoint->GetPhysicalVolume();
159  std::string nameVolume = ForwardName::getName(currentPV->GetName());
160  std::string postnameVolume = ForwardName::getName(postPV->GetName());
161  G4String particleType = theTrack->GetDefinition()->GetParticleName();
162  edm::LogVerbatim("ForwardSim") << "ZdcSD:: getEnergyDeposit: \n"
163  << " preStepPoint: " << nameVolume << "," << stepL << "," << stepE << "," << beta
164  << "," << charge << "\n"
165  << " postStepPoint: " << postnameVolume << "," << costheta << "," << theta << ","
166  << eta << "," << phi << "," << particleType << " id= " << theTrack->GetTrackID()
167  << " Etot(GeV)= " << theTrack->GetTotalEnergy() / GeV;
168 #endif
169  const double bThreshold = 0.67;
170  if (beta > bThreshold) {
171 #ifdef EDM_ML_DEBUG
172  edm::LogVerbatim("ForwardSim") << "ZdcSD:: getEnergyDeposit: pass ";
173 #endif
174  const float nMedium = 1.4925;
175  // float photEnSpectrDL = 10714.285714;
176  // photEnSpectrDL = (1./400.nm-1./700.nm)*10000000.cm/nm; /* cm-1 */
177 
178  const float photEnSpectrDE = 1.24;
179  // E = 2pi*(1./137.)*(eV*cm/370.)/lambda = 12.389184*(eV*cm)/lambda
180  // Emax = 12.389184*(eV*cm)/400nm*10-7cm/nm = 3.01 eV
181  // Emin = 12.389184*(eV*cm)/700nm*10-7cm/nm = 1.77 eV
182  // delE = Emax - Emin = 1.24 eV
183 
184  const float effPMTandTransport = 0.15;
185 
186  // Check these values
187  const float thFullRefl = 23.;
188  float thFullReflRad = thFullRefl * pi / 180.;
189 
190  float thFibDirRad = thFibDir * pi / 180.;
191 
192  // at which theta the point is located:
193  // float th1 = hitPoint.theta();
194 
195  // theta of charged particle in LabRF(hit momentum direction):
196  float costh = hit_mom.z() / sqrt(hit_mom.x() * hit_mom.x() + hit_mom.y() * hit_mom.y() + hit_mom.z() * hit_mom.z());
197  float th = acos(std::min(std::max(costh, -1.f), 1.f));
198  // just in case (can do both standard ranges of phi):
199  if (th < 0.)
200  th += twopi;
201 
202  // theta of cone with Cherenkov photons w.r.t.direction of charged part.:
203  float costhcher = 1. / (nMedium * beta);
204  float thcher = acos(std::min(std::max(costhcher, -1.f), 1.f));
205 
206  // diff thetas of charged part. and quartz direction in LabRF:
207  float DelFibPart = std::abs(th - thFibDirRad);
208 
209  // define real distances:
210  float d = std::abs(std::tan(th) - std::tan(thFibDirRad));
211 
212  float a = std::tan(thFibDirRad) + std::tan(std::abs(thFibDirRad - thFullReflRad));
213  float r = std::tan(th) + std::tan(std::abs(th - thcher));
214 
215  // define losses d_qz in cone of full reflection inside quartz direction
216  float d_qz = -1;
217 #ifdef EDM_ML_DEBUG
218  float variant = -1;
219 #endif
220  // if (d > (r+a))
221  if (DelFibPart > (thFullReflRad + thcher)) {
222 #ifdef EDM_ML_DEBUG
223  variant = 0.;
224 #endif
225  d_qz = 0.;
226  } else {
227  // if ((DelFibPart + thcher) < thFullReflRad ) [(d+r) < a]
228  if ((th + thcher) < (thFibDirRad + thFullReflRad) && (th - thcher) > (thFibDirRad - thFullReflRad)) {
229 #ifdef EDM_ML_DEBUG
230  variant = 1.;
231 #endif
232  d_qz = 1.;
233  } else {
234  // if ((thcher - DelFibPart ) > thFullReflRad ) [(r-d) > a]
235  if ((thFibDirRad + thFullReflRad) < (th + thcher) && (thFibDirRad - thFullReflRad) > (th - thcher)) {
236 #ifdef EDM_ML_DEBUG
237  variant = 2.;
238 #endif
239  d_qz = 0.;
240  } else {
241 #ifdef EDM_ML_DEBUG
242  variant = 3.; // d_qz is calculated below
243 #endif
244  // use crossed length of circles(cone projection) - dC1/dC2 :
245  float arg_arcos = 0.;
246  float tan_arcos = 2. * a * d;
247  if (tan_arcos != 0.)
248  arg_arcos = (r * r - a * a - d * d) / tan_arcos;
249  arg_arcos = std::abs(arg_arcos);
250  float th_arcos = acos(std::min(std::max(arg_arcos, -1.f), 1.f));
251  d_qz = th_arcos / twopi;
252  d_qz = std::abs(d_qz);
253 #ifdef EDM_ML_DEBUG
254  edm::LogVerbatim("ForwardSim") << " d_qz: " << r << "," << a << "," << d << " " << tan_arcos << " "
255  << arg_arcos;
256  edm::LogVerbatim("ForwardSim") << "," << arg_arcos;
257  edm::LogVerbatim("ForwardSim") << " " << d_qz;
258  edm::LogVerbatim("ForwardSim") << " " << th_arcos;
259  edm::LogVerbatim("ForwardSim") << "," << d_qz;
260 #endif
261  }
262  }
263  }
264  double meanNCherPhot = 0.;
265  int poissNCherPhot = 0;
266  if (d_qz > 0) {
267  meanNCherPhot = 370. * charge * charge * (1. - 1. / (nMedium * nMedium * beta * beta)) * photEnSpectrDE * stepL;
268 
269  poissNCherPhot = std::max((int)G4Poisson(meanNCherPhot), 0);
270  NCherPhot = poissNCherPhot * effPMTandTransport * d_qz;
271  }
272 
273 #ifdef EDM_ML_DEBUG
274  edm::LogVerbatim("ForwardSim") << "ZdcSD:: getEnergyDeposit: gED: " << stepE << "," << costh << "," << th << ","
275  << costhcher << "," << thcher << "," << DelFibPart << "," << d << "," << a << ","
276  << r << "," << hitPoint << "," << hit_mom << "," << vert_mom << "," << localPoint
277  << "," << charge << "," << beta << "," << stepL << "," << d_qz << "," << variant
278  << "," << meanNCherPhot << "," << poissNCherPhot << "," << NCherPhot;
279 #endif
280  // --constants-----------------
281  // << "," << photEnSpectrDE
282  // << "," << nMedium
283  // << "," << bThreshold
284  // << "," << thFibDirRad
285  // << "," << thFullReflRad
286  // << "," << effPMTandTransport
287  // --other variables-----------
288  // << "," << curprocess
289  // << "," << nameProcess
290  // << "," << name
291  // << "," << rad
292  // << "," << mat
293 
294  } else {
295  // determine failure mode: beta, charge, and/or nameVolume
296  if (beta <= bThreshold)
297  edm::LogVerbatim("ForwardSim") << "ZdcSD:: getEnergyDeposit: fail beta=" << beta;
298 
299 #ifdef EDM_ML_DEBUG
300  std::string nameVolume = ForwardName::getName(currentPV->GetName());
301  if (!(nameVolume == "ZDC_EMFiber" || nameVolume == "ZDC_HadFiber"))
302  edm::LogVerbatim("ForwardSim") << "ZdcSD:: getEnergyDeposit: fail nv=" << nameVolume;
303 #endif
304  }
305 
306  return NCherPhot;
307 }
308 
309 uint32_t ZdcSD::setDetUnitId(const G4Step* aStep) {
310  return (numberingScheme.get() == nullptr ? 0 : numberingScheme.get()->getUnitID(aStep));
311 }
312 
314  if (scheme != nullptr) {
315  edm::LogVerbatim("ForwardSim") << "ZdcSD: updates numbering scheme for " << GetName();
316  numberingScheme.reset(scheme);
317  }
318 }
float edepositEM
Definition: CaloSD.h:140
void setNumberingScheme(ZdcNumberingScheme *scheme)
Definition: ZdcSD.cc:313
Log< level::Info, true > LogVerbatim
double thFibDir
Definition: ZdcSD.h:33
int getTrackID() const
Definition: CaloG4Hit.h:64
std::unique_ptr< ZdcNumberingScheme > numberingScheme
Definition: ZdcSD.h:37
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
math::XYZPoint getEntryLocal() const
Definition: CaloG4Hit.h:49
Definition: CaloSD.h:40
int verbosity
Definition: ZdcSD.h:31
double getEM() const
Definition: CaloG4Hit.h:55
std::vector< ZdcShowerLibrary::Hit > hits
Definition: ZdcSD.h:38
G4ThreeVector posGlobal
Definition: CaloSD.h:138
bool useShowerHits
Definition: ZdcSD.h:32
void processHit(const G4Step *step)
Definition: CaloSD.h:113
double getHadr() const
Definition: CaloG4Hit.h:58
const Double_t pi
float edepositHAD
Definition: CaloSD.h:140
T sqrt(T t)
Definition: SSEVec.h:19
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:653
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
CaloG4Hit * currentHit
Definition: CaloSD.h:146
d
Definition: ztail.py:151
ZdcSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: ZdcSD.cc:30
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID
Definition: CaloSD.h:142
void initRun() override
Definition: ZdcSD.cc:64
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
std::unique_ptr< ZdcShowerLibrary > showerLibrary
Definition: ZdcSD.h:36
float incidentEnergy
Definition: CaloSD.h:139
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:826
double a
Definition: hdecay.h:119
double zdcHitEnergyCut
Definition: ZdcSD.h:34
std::string getName(const G4String &)
Definition: ForwardName.cc:3
double getIncidentEnergy() const
Definition: CaloG4Hit.h:61
bool getFromLibrary(const G4Step *) override
Definition: ZdcSD.cc:66
Geom::Theta< T > theta() const
uint32_t setDetUnitId(const G4Step *step) override
Definition: ZdcSD.cc:309
bool useShowerLibrary
Definition: ZdcSD.h:32
G4ThreeVector entrancePoint
Definition: CaloSD.h:136
G4ThreeVector entranceLocal
Definition: CaloSD.h:137
double getEnergyDeposit(const G4Step *) override
Definition: ZdcSD.cc:120
void setParameterized(bool val)
Definition: CaloSD.h:110