CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
ZdcSD Class Reference

#include <ZdcSD.h>

Inheritance diagram for ZdcSD:
CaloSD SensitiveCaloDetector Observer< const BeginOfRun *> Observer< const BeginOfEvent *> Observer< const BeginOfTrack *> Observer< const EndOfTrack *> Observer< const EndOfEvent *> SensitiveDetector

Public Member Functions

double calculateCherenkovDeposit (const G4Step *)
 
double calculateMeanNumberOfPhotons (int, double, double)
 
double calculateN2InvIntegral (double)
 
double convertEnergyToWavelength (double)
 
double evaluateFunction (const std::vector< double > &, const std::vector< double > &, double)
 
double generatePhotonEnergy (int, double, double)
 
double linearInterpolation (double, double, double, double, double)
 
double photonEnergyDist (int, double, double)
 
double pmtEfficiency (double)
 
bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory) override
 
uint32_t setDetUnitId (const G4Step *step) override
 
 ZdcSD (const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
 ~ZdcSD () override=default
 
- Public Member Functions inherited from CaloSD
 CaloSD (const std::string &aSDname, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, float timeSlice=1., bool ignoreTkID=false, const std::string &newcolname="")
 
void clear () override
 
void clearHits () override
 
void DrawAll () override
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *HCE) override
 
bool isItFineCalo (const G4VTouchable *touch)
 
void newCollection (const std::string &name, edm::ParameterSet const &p)
 
void PrintAll () override
 
bool ProcessHits (G4GFlashSpot *aSpot, G4TouchableHistory *) override
 
void reset () override
 
 ~CaloSD () override
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (const std::string &iname, const SensitiveDetectorCatalog &clg, const std::string &newcollname="")
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const SensitiveDetectorCatalog &, bool calo, const std::string &newcollname="")
 
 ~SensitiveDetector () override
 
- Public Member Functions inherited from Observer< const BeginOfRun *>
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent *>
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack *>
 Observer ()
 
void slotForUpdate (const BeginOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack *>
 Observer ()
 
void slotForUpdate (const EndOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent *>
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 

Protected Member Functions

void initRun () override
 
- Protected Member Functions inherited from CaloSD
bool checkHit (int k=0)
 
CaloG4HitcreateNewHit (const G4Step *, const G4Track *, int k)
 
virtual void endEvent ()
 
virtual double EnergyCorrected (const G4Step &step, const G4Track *)
 
virtual bool filterHit (CaloG4Hit *, double)
 
unsigned int findBoundaryCrossingParent (const G4Track *track, bool markParentAsSaveable=true)
 
double getAttenuation (const G4Step *aStep, double birk1, double birk2, double birk3) const
 
virtual uint16_t getDepth (const G4Step *)
 
virtual double getEnergyDeposit (const G4Step *step)
 
virtual bool getFromLibrary (const G4Step *step)
 
int getNumberOfHits (int k=0)
 
double getResponseWt (const G4Track *, int k=0)
 
virtual int getTrackID (const G4Track *)
 
bool hitExists (const G4Step *, int k)
 
void ignoreRejection ()
 
virtual void initEvent (const BeginOfEvent *)
 
void printDetectorLevels (const G4VTouchable *) const
 
void processHit (const G4Step *step)
 
virtual void processSecondHit (const G4Step *, const G4Track *)
 
void resetForNewPrimary (const G4Step *)
 
void setNumberCheckedHits (int val, int k=0)
 
void setParameterized (bool val)
 
G4ThreeVector setToGlobal (const G4ThreeVector &, const G4VTouchable *) const
 
G4ThreeVector setToLocal (const G4ThreeVector &, const G4VTouchable *) const
 
void setUseMap (bool val)
 
std::string shortreprID (const CaloHitID &ID)
 
std::string shortreprID (const CaloG4Hit *hit)
 
void update (const BeginOfRun *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfEvent *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const ::EndOfEvent *) override
 
void updateHit (CaloG4Hit *, int k)
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 
- Protected Member Functions inherited from Observer< const EndOfEvent *>
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 

Private Member Functions

int setTrackID (const G4Step *step) override
 

Private Attributes

std::vector< ZdcShowerLibrary::Hithits
 
std::unique_ptr< ZdcNumberingSchemenumberingScheme
 
std::unique_ptr< ZdcShowerLibraryshowerLibrary
 
double thFibDir
 
bool useShowerHits
 
bool useShowerLibrary
 
int verbosity
 
double zdcHitEnergyCut
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Static Protected Member Functions inherited from CaloSD
static std::string printableDecayChain (const std::vector< unsigned int > &decayChain)
 
- Protected Attributes inherited from CaloSD
std::string collName_ [2]
 
CaloG4HitcurrentHit [2]
 
CaloHitID currentID [2]
 
std::string detName_
 
float edepositEM
 
float edepositHAD
 
double eminHit
 
double energyCut
 
G4ThreeVector entranceLocal
 
G4ThreeVector entrancePoint
 
bool forceSave
 
std::vector< std::string > hcn_
 
float incidentEnergy
 
double kmaxIon
 
double kmaxNeutron
 
double kmaxProton
 
int nHC_
 
G4ThreeVector posGlobal
 
CaloHitID previousID [2]
 
bool suppressHeavy
 
double tmaxHit
 
std::vector< int > useResMap_
 

Detailed Description

Definition at line 13 of file ZdcSD.h.

Constructor & Destructor Documentation

◆ ZdcSD()

ZdcSD::ZdcSD ( const std::string &  name,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 32 of file ZdcSD.cc.

References edm::ParameterSet::getParameter(), Skims_PA_cff::name, numberingScheme, AlCaHLTBitMon_ParallelJobs::p, CaloSD::setParameterized(), showerLibrary, useShowerHits, useShowerLibrary, verbosity, and zdcHitEnergyCut.

36  : CaloSD(name, clg, p, manager) {
37  edm::ParameterSet m_ZdcSD = p.getParameter<edm::ParameterSet>("ZdcSD");
38  useShowerLibrary = m_ZdcSD.getParameter<bool>("UseShowerLibrary");
39  useShowerHits = m_ZdcSD.getParameter<bool>("UseShowerHits");
40  zdcHitEnergyCut = m_ZdcSD.getParameter<double>("ZdcHitEnergyCut") * GeV;
41  verbosity = m_ZdcSD.getParameter<int>("Verbosity");
42  int verbn = verbosity / 10;
43  verbosity %= 10;
44  numberingScheme = std::make_unique<ZdcNumberingScheme>(verbn);
45 
46  edm::LogVerbatim("ForwardSim") << "***************************************************\n"
47  << "* *\n"
48  << "* Constructing a ZdcSD with name " << name << " *\n"
49  << "* *\n"
50  << "***************************************************";
51 
52  edm::LogVerbatim("ForwardSim") << "\nUse of shower library is set to " << useShowerLibrary
53  << "\nUse of Shower hits method is set to " << useShowerHits;
54 
55  edm::LogVerbatim("ForwardSim") << "\nEnergy Threshold Cut set to " << zdcHitEnergyCut / CLHEP::GeV << " (GeV)";
56 
57  if (useShowerLibrary) {
58  showerLibrary = std::make_unique<ZdcShowerLibrary>(name, p);
59  setParameterized(true);
60  } else {
61  showerLibrary.reset(nullptr);
62  }
63 }
Log< level::Info, true > LogVerbatim
std::unique_ptr< ZdcNumberingScheme > numberingScheme
Definition: ZdcSD.h:44
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
int verbosity
Definition: ZdcSD.h:38
bool useShowerHits
Definition: ZdcSD.h:39
CaloSD(const std::string &aSDname, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, float timeSlice=1., bool ignoreTkID=false, const std::string &newcolname="")
Definition: CaloSD.cc:33
std::unique_ptr< ZdcShowerLibrary > showerLibrary
Definition: ZdcSD.h:43
double zdcHitEnergyCut
Definition: ZdcSD.h:42
bool useShowerLibrary
Definition: ZdcSD.h:39
void setParameterized(bool val)
Definition: CaloSD.h:114

◆ ~ZdcSD()

ZdcSD::~ZdcSD ( )
overridedefault

Member Function Documentation

◆ calculateCherenkovDeposit()

double ZdcSD::calculateCherenkovDeposit ( const G4Step *  aStep)

Definition at line 131 of file ZdcSD.cc.

References a, b, HLT_2023v12_cff::beta, HltBtagPostValidation_cff::c, calculateMeanNumberOfPhotons(), ALCARECOTkAlJpsiMuMu_cff::charge, convertEnergyToWavelength(), funct::cos(), EMAX, EMIN, mps_fire::i, mag(), dqmiodumpmetadata::n, NAperRINDEX, run3scouting_cff::nPhotons, pmtEfficiency(), findAndChange::post, findAndChange::pre, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, dttmaxenums::R, RINDEX, CaloSD::setToLocal(), funct::sin(), mathSSE::sqrt(), submitPVValidationJobs::t, and V0Monitor_cff::v0.

Referenced by ProcessHits().

131  {
132  G4Material* material = aStep->GetTrack()->GetMaterial();
133 
134  if (material->GetName() != "quartz")
135  return 0.0; // 0 deposit if material is not quartz
136  else {
137  const G4StepPoint* pPreStepPoint = aStep->GetPreStepPoint();
138  const G4StepPoint* pPostStepPoint = aStep->GetPostStepPoint();
139  const G4String volumeName = pPreStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume()->GetName();
140 
141  G4ThreeVector pre = pPreStepPoint->GetPosition();
142  G4ThreeVector post = pPostStepPoint->GetPosition();
143 
144  if ((post - pre).mag() < 1E-9)
145  return 0.0;
146 
147  //Convert step coordinates to local (fiber) coodinates
148  const G4ThreeVector localPre = setToLocal(pre, pPreStepPoint->GetTouchable());
149  const G4ThreeVector localPost = setToLocal(post, pPreStepPoint->GetTouchable());
150  // Calculate the unit direction vector in local coordinates
151 
152  const G4ThreeVector particleDirection = (localPost - localPre) / (localPost - localPre).mag();
153 
154  const G4DynamicParticle* aParticle = aStep->GetTrack()->GetDynamicParticle();
155  int charge = round(aParticle->GetDefinition()->GetPDGCharge());
156 
157  if (charge == 0)
158  return 0.0;
159 
160  double beta = 0.5 * (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta());
161  double stepLength = aStep->GetStepLength() / 1000; // Geant4 stepLength is in "mm"
162 
163  int nPhotons; // Number of Cherenkov photons
164 
165  nPhotons = G4Poisson(calculateMeanNumberOfPhotons(charge, beta, stepLength));
166 
167  double totalE = 0.0;
168 
169  for (int i = 0; i < nPhotons; i++) {
170  // uniform refractive index in PMT range -> uniform energy distribution
171  double photonE = EMIN + G4UniformRand() * (EMAX - EMIN);
172  // UPDATE: taking into account dispersion relation -> energy distribution
173 
174  if (G4UniformRand() > pmtEfficiency(convertEnergyToWavelength(photonE)))
175  continue;
176 
177  double omega = G4UniformRand() * twopi;
178  double thetaC = acos(1.0 / (beta * RINDEX));
179 
180 #ifdef EDM_ML_DEBUG
181  edm::LogVerbatim("ZdcSD") << "E_gamma: " << photonE << "\t omega: " << omega << "\t thetaC: " << thetaC;
182 #endif
183  // Calculate momentum direction w.r.t primary particle (z-direction)
184  double px = photonE * sin(thetaC) * cos(omega);
185  double py = photonE * sin(thetaC) * sin(omega);
186  double pz = photonE * cos(thetaC);
187  G4ThreeVector photonMomentum(px, py, pz);
188 
189 #ifdef EDM_ML_DEBUG
190  edm::LogVerbatim("ZdcSD") << "pPR = (" << particleDirection.x() << "," << particleDirection.y() << ","
191  << particleDirection.z() << ")";
192  edm::LogVerbatim("ZdcSD") << "pCH = (" << px << "," << py << "," << pz << ")";
193 #endif
194  // Rotate to the fiber reference frame
195  photonMomentum.rotateUz(particleDirection);
196 
197 #ifdef EDM_ML_DEBUG
198  edm::LogVerbatim("ZdcSD") << "pLAB = (" << photonMomentum.x() << "," << photonMomentum.y() << ","
199  << photonMomentum.z() << ")";
200 #endif
201  // Get random position along G4Step
202  G4ThreeVector photonPosition = localPre + G4UniformRand() * (localPost - localPre);
203 
204  // 2D vectors to calculate impact position (x*,y*)
205  G4TwoVector r0(photonPosition);
206  G4TwoVector v0(photonMomentum);
207 
208  double R = 0.3; /*mm, fiber radius*/
209  double R2 = 0.3 * 0.3;
210 
211  if (r0.mag() < R && photonMomentum.z() < 0.0) {
212  // 2nd order polynomial coefficients
213  double a = v0.mag2();
214  double b = 2.0 * r0 * v0;
215  double c = r0.mag2() - R2;
216 
217  if (a < 1E-6)
218  totalE += 1; //photonE /*eV*/;
219  else {
220  // calculate intersection point - solving 2nd order polynomial
221  double t = (-b + sqrt(b * b - 4.0 * a * c)) / (2.0 * a);
222  G4ThreeVector n(r0.x() + v0.x() * t, r0.y() + v0.y() * t, 0.0); // surface normal
223  double cosTheta = (n * photonMomentum) / (n.mag() * photonE); // cosine of incident angle
224 
225  if (cosTheta >= NAperRINDEX) // lightguide condition
226  totalE += 1; //photonE /*eV*/;
227  }
228  }
229 
230 #ifdef EDM_ML_DEBUG
231  edm::LogVerbatim("ZdcSD") << "r = (" << photonPosition.x() << "," << photonPosition.y() << ","
232  << photonPosition.z() << ")" << std::endl;
233 #endif
234  }
235 
236 #ifdef EDM_ML_DEBUG
237  if (nPhotons > 30) {
238  edm::LogVerbatim("ZdcSD") << totalE;
239 
240  if (totalE > 0)
241  edm::LogVerbatim("ZdcSD") << pre.x() << " " << pre.y() << " " << pre.z() << " " << totalE << std::endl;
242  }
243 #endif
244  return totalE;
245  }
246 }
double convertEnergyToWavelength(double)
Definition: ZdcSD.cc:357
Log< level::Info, true > LogVerbatim
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
Definition: CaloSD.cc:454
const double EMIN
Definition: ZdcSD.cc:126
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
double b
Definition: hdecay.h:120
double calculateMeanNumberOfPhotons(int, double, double)
Definition: ZdcSD.cc:250
double pmtEfficiency(double)
Definition: ZdcSD.cc:299
double a
Definition: hdecay.h:121
const double RINDEX
Definition: ZdcSD.cc:122
const double NAperRINDEX
Definition: ZdcSD.cc:124
const double EMAX
Definition: ZdcSD.cc:125

◆ calculateMeanNumberOfPhotons()

double ZdcSD::calculateMeanNumberOfPhotons ( int  charge,
double  beta,
double  stepLength 
)

Definition at line 250 of file ZdcSD.cc.

References ALPHA, HLT_2023v12_cff::beta, ALCARECOTkAlJpsiMuMu_cff::charge, EMAX, EMIN, HBARC, and RINDEX.

Referenced by calculateCherenkovDeposit().

250  {
251  // Return mean number of Cherenkov photons
252  return (ALPHA * charge * charge * stepLength) / HBARC * (EMAX - EMIN) * (1.0 - 1.0 / (beta * beta * RINDEX * RINDEX));
253 }
const double HBARC
Definition: ZdcSD.cc:128
const double EMIN
Definition: ZdcSD.cc:126
const double ALPHA
Definition: ZdcSD.cc:127
const double RINDEX
Definition: ZdcSD.cc:122
const double EMAX
Definition: ZdcSD.cc:125

◆ calculateN2InvIntegral()

double ZdcSD::calculateN2InvIntegral ( double  Emin)

Definition at line 284 of file ZdcSD.cc.

References evaluateFunction().

284  {
285  // Hardcoded minimum photon energy table (eV)
286  const std::vector<double> EMIN_TAB{1.75715, 1.81902, 1.88311, 1.94944, 2.0183, 2.08939, 2.16302, 2.23919,
287  2.31789, 2.39954, 2.48416, 2.57175, 2.66232, 2.75643, 2.85349, 2.95411,
288  3.05756, 3.16528, 3.2774, 3.39218, 3.5123, 3.6359, 3.76394, 3.89642,
289  4.03332, 4.17596, 4.32302, 4.47596, 4.63319};
290 
291  // Hardcoded integral values
292  const std::vector<double> INTEGRAL_TAB{1.39756, 1.36835, 1.33812, 1.30686, 1.27443, 1.24099, 1.20638, 1.17061,
293  1.1337, 1.09545, 1.05586, 1.01493, 0.972664, 0.928815, 0.883664, 0.836938,
294  0.788988, 0.739157, 0.687404, 0.634547, 0.579368, 0.522743, 0.464256, 0.40393,
295  0.341808, 0.27732, 0.211101, 0.142536, 0.0723891};
296  return evaluateFunction(EMIN_TAB, INTEGRAL_TAB, Emin);
297 }
double evaluateFunction(const std::vector< double > &, const std::vector< double > &, double)
Definition: ZdcSD.cc:329

◆ convertEnergyToWavelength()

double ZdcSD::convertEnergyToWavelength ( double  energy)

Definition at line 357 of file ZdcSD.cc.

References hcalRecHitTable_cff::energy.

Referenced by calculateCherenkovDeposit().

357 { return (1240.0 / energy); }

◆ evaluateFunction()

double ZdcSD::evaluateFunction ( const std::vector< double > &  X,
const std::vector< double > &  Y,
double  x 
)

Definition at line 329 of file ZdcSD.cc.

References mps_fire::i, linearInterpolation(), x, X, and beamSpotPI::Y.

Referenced by calculateN2InvIntegral(), and photonEnergyDist().

329  {
330  for (unsigned int i = 0; i < X.size() - 1; i++) {
331  if (x > X[i] && x < X[i + 1]) {
332 #ifdef EDM_ML_DEBUG
333  edm::LogVerbatim("ZdcSD") << X[i] << "\t" << Y[i] << "\t" << X[i + 1] << "\t" << Y[i + 1] << "\t" << x << "\t"
334  << linearInterpolation(X[i], Y[i], X[i + 1], Y[i + 1], x);
335 #endif
336  return linearInterpolation(X[i], Y[i], X[i + 1], Y[i + 1], x);
337  }
338  }
339 
340  if (fabs(X[0] - x) < 1E-8)
341  return Y[0];
342  else if (fabs(X[X.size() - 1] - x) < 1E-8)
343  return Y[X.size() - 1];
344  else
345  return NAN;
346 }
Log< level::Info, true > LogVerbatim
#define X(str)
Definition: MuonsGrabber.cc:38
double linearInterpolation(double, double, double, double, double)
Definition: ZdcSD.cc:349

◆ generatePhotonEnergy()

double ZdcSD::generatePhotonEnergy ( int  charge,
double  beta,
double  Emin 
)

Definition at line 271 of file ZdcSD.cc.

References HLT_2023v12_cff::beta, ALCARECOTkAlJpsiMuMu_cff::charge, EMAX, and photonEnergyDist().

271  {
272  double photonE;
273 
274  // Use rejection method
275  do {
276  photonE = G4UniformRand() * (EMAX - Emin) + Emin;
277  } while (G4UniformRand() > photonEnergyDist(photonE, charge, beta) / photonEnergyDist(EMAX, charge, beta));
278 
279  return photonE;
280 }
double photonEnergyDist(int, double, double)
Definition: ZdcSD.cc:256
const double EMAX
Definition: ZdcSD.cc:125

◆ initRun()

void ZdcSD::initRun ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 65 of file ZdcSD.cc.

References hits, showerLibrary, and useShowerLibrary.

65  {
66  if (useShowerLibrary) {
67  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
68  showerLibrary->initRun(theParticleTable);
69  }
70  hits.clear();
71 }
std::vector< ZdcShowerLibrary::Hit > hits
Definition: ZdcSD.h:46
std::unique_ptr< ZdcShowerLibrary > showerLibrary
Definition: ZdcSD.h:43
bool useShowerLibrary
Definition: ZdcSD.h:39

◆ linearInterpolation()

double ZdcSD::linearInterpolation ( double  x1,
double  y1,
double  x2,
double  y2,
double  z 
)

◆ photonEnergyDist()

double ZdcSD::photonEnergyDist ( int  charge,
double  beta,
double  E 
)

Definition at line 256 of file ZdcSD.cc.

References ALPHA, HLT_2023v12_cff::beta, ALCARECOTkAlJpsiMuMu_cff::charge, evaluateFunction(), and HBARC.

Referenced by generatePhotonEnergy().

256  {
257  const std::vector<double> ENERGY_TAB{1.75715, 1.81902, 1.88311, 1.94944, 2.0183, 2.08939, 2.16302, 2.23919,
258  2.31789, 2.39954, 2.48416, 2.57175, 2.66232, 2.75643, 2.85349, 2.95411,
259  3.05756, 3.16528, 3.2774, 3.39218, 3.5123, 3.6359, 3.76394, 3.89642,
260  4.03332, 4.17596, 4.32302, 4.47596, 4.63319, 4.79629};
261 
262  const std::vector<double> RINDEX_TAB{1.45517, 1.45572, 1.45631, 1.45693, 1.45758, 1.45826, 1.45899, 1.45976,
263  1.46057, 1.46144, 1.46238, 1.46337, 1.46444, 1.46558, 1.4668, 1.46812,
264  1.46953, 1.47105, 1.4727, 1.47447, 1.4764, 1.47847, 1.48071, 1.48315,
265  1.48579, 1.48868, 1.49182, 1.49526, 1.499, 1.5031};
266  double rIndex = evaluateFunction(ENERGY_TAB, RINDEX_TAB, E);
267  return (ALPHA * charge * charge) / HBARC * (1.0 - 1.0 / (beta * beta * rIndex * rIndex));
268 }
const double HBARC
Definition: ZdcSD.cc:128
double evaluateFunction(const std::vector< double > &, const std::vector< double > &, double)
Definition: ZdcSD.cc:329
const double ALPHA
Definition: ZdcSD.cc:127

◆ pmtEfficiency()

double ZdcSD::pmtEfficiency ( double  lambda)

Definition at line 299 of file ZdcSD.cc.

References a, b, and mps_fire::i.

Referenced by calculateCherenkovDeposit().

299  {
300  // Hardcoded wavelength values (nm)
301  const std::vector<double> LAMBDA_TAB{263.27, 265.98, 268.69, 271.39, 273.20, 275.90, 282.22, 282.22, 293.04,
302  308.38, 325.52, 346.26, 367.91, 392.27, 417.53, 440.98, 463.53, 484.28,
303  502.32, 516.75, 528.48, 539.30, 551.93, 564.56, 574.48, 584.41, 595.23,
304  606.96, 616.88, 625.00, 632.22, 637.63, 642.14, 647.55, 652.96, 656.57,
305  661.08, 666.49, 669.20, 673.71, 677.32, 680.93, 686.34, 692.65};
306 
307  // Hardcoded quantum efficiency values
308  const std::vector<double> EFF_TAB{2.215, 2.860, 3.659, 4.724, 5.989, 7.734, 9.806, 9.806, 12.322,
309  15.068, 17.929, 20.570, 22.963, 24.050, 23.847, 22.798, 20.445, 18.003,
310  15.007, 12.282, 9.869, 7.858, 6.373, 5.121, 4.077, 3.276, 2.562,
311  2.077, 1.669, 1.305, 1.030, 0.805, 0.629, 0.492, 0.388, 0.303,
312  0.239, 0.187, 0.144, 0.113, 0.088, 0.069, 0.054, 0.043};
313  //double efficiency = evaluateFunction(LAMBDA_TAB,EFF_TAB,lambda);
314 
315  // Using linear interpolation to calculate efficiency
316  for (int i = 0; i < 44 - 1; i++) {
317  if (lambda > LAMBDA_TAB[i] && lambda < LAMBDA_TAB[i + 1]) {
318  double a = (EFF_TAB[i] - EFF_TAB[i + 1]) / (LAMBDA_TAB[i] - LAMBDA_TAB[i + 1]);
319  double b = EFF_TAB[i] - a * LAMBDA_TAB[i];
320  return (a * lambda + b) / 100.0;
321  }
322  }
323 
324  return 0;
325 }
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121

◆ ProcessHits()

bool ZdcSD::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
)
overridevirtual

Reimplemented from CaloSD.

Definition at line 73 of file ZdcSD.cc.

References calculateCherenkovDeposit(), CaloSD::createNewHit(), CaloSD::currentHit, CaloSD::currentID, hcalRecHitTable_cff::depth, CaloSD::edepositEM, CaloSD::edepositHAD, hcalRecHitTable_cff::energy, CaloSD::getDepth(), CaloSD::getTrackID(), CaloSD::hitExists(), G4TrackToParticleID::isGammaElectronPositron(), SensitiveDetector::NaNTrap(), findAndChange::pre, setDetUnitId(), CaloHitID::setID(), hcalRecHitTable_cff::time, and useShowerHits.

73  {
74  NaNTrap(aStep);
75 
76  /*
77  if (useShowerLibrary)
78  getFromLibrary(aStep);
79  */
80 #ifdef EDM_ML_DEBUG
81  edm::LogVerbatim("ZdcSD") << "ZdcSD::" << GetName() << " ID= " << aStep->GetTrack()->GetTrackID()
82  << " prID= " << aStep->GetTrack()->GetParentID()
83  << " Eprestep= " << aStep->GetPreStepPoint()->GetKineticEnergy()
84  << " step= " << aStep->GetStepLength() << " Edep= " << aStep->GetTotalEnergyDeposit();
85 #endif
86  if (useShowerHits) {
87  // check unitID
88  unsigned int unitID = setDetUnitId(aStep);
89  if (unitID == 0)
90  return false;
91 
92  auto const theTrack = aStep->GetTrack();
93  uint16_t depth = getDepth(aStep);
94 
95  double time = theTrack->GetGlobalTime() / nanosecond;
96  int primaryID = getTrackID(theTrack);
97  currentID[0].setID(unitID, time, primaryID, depth);
98  double energy = calculateCherenkovDeposit(aStep);
100  edepositEM = energy;
101  edepositHAD = 0;
102  } else {
103  edepositEM = 0;
105  }
106  if (!hitExists(aStep, 0) && edepositEM + edepositHAD > 0.) {
107 #ifdef EDM_ML_DEBUG
108  G4ThreeVector pre = aStep->GetPreStepPoint()->GetPosition();
109  edm::LogVerbatim("ZdcSD") << pre.x() << " " << pre.y() << " " << pre.z();
110 #endif
111  currentHit[0] = CaloSD::createNewHit(aStep, theTrack, 0);
112  }
113  }
114  return true;
115 }
CaloG4Hit * currentHit[2]
Definition: CaloSD.h:152
float edepositEM
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
virtual uint16_t getDepth(const G4Step *)
Definition: CaloSD.cc:912
bool useShowerHits
Definition: ZdcSD.h:39
virtual int getTrackID(const G4Track *)
Definition: CaloSD.cc:872
double calculateCherenkovDeposit(const G4Step *)
Definition: ZdcSD.cc:131
CaloG4Hit * createNewHit(const G4Step *, const G4Track *, int k)
Definition: CaloSD.cc:623
float edepositHAD
Definition: CaloSD.h:144
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
bool hitExists(const G4Step *, int k)
Definition: CaloSD.cc:462
static bool isGammaElectronPositron(int pdgCode)
CaloHitID currentID[2]
Definition: CaloSD.h:146
uint32_t setDetUnitId(const G4Step *step) override
Definition: ZdcSD.cc:361
void NaNTrap(const G4Step *step) const

◆ setDetUnitId()

uint32_t ZdcSD::setDetUnitId ( const G4Step *  step)
overridevirtual

Implements CaloSD.

Definition at line 361 of file ZdcSD.cc.

References numberingScheme.

Referenced by ProcessHits().

361  {
362  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
363 }
std::unique_ptr< ZdcNumberingScheme > numberingScheme
Definition: ZdcSD.h:44

◆ setTrackID()

int ZdcSD::setTrackID ( const G4Step *  step)
overrideprivatevirtual

Reimplemented from CaloSD.

Definition at line 365 of file ZdcSD.cc.

References TrackInformation::getIDonCaloSurface(), CaloSD::previousID, and CaloSD::resetForNewPrimary().

365  {
366  auto const theTrack = aStep->GetTrack();
367  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
368  int primaryID = trkInfo->getIDonCaloSurface();
369  if (primaryID == 0) {
370 #ifdef EDM_ML_DEBUG
371  auto const preStepPoint = aStep->GetPreStepPoint();
372  double etrack = preStepPoint->GetKineticEnergy();
373  edm::LogVerbatim("ZdcSD") << "ZdcSD: Problem with primaryID **** set by force to TkID **** "
374  << theTrack->GetTrackID() << " E " << etrack;
375 #endif
376  primaryID = theTrack->GetTrackID();
377  }
378  if (primaryID != previousID[0].trackID())
379  resetForNewPrimary(aStep);
380  return primaryID;
381 }
Log< level::Info, true > LogVerbatim
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:713
int getIDonCaloSurface() const
CaloHitID previousID[2]
Definition: CaloSD.h:146

Member Data Documentation

◆ hits

std::vector<ZdcShowerLibrary::Hit> ZdcSD::hits
private

Definition at line 46 of file ZdcSD.h.

Referenced by initRun().

◆ numberingScheme

std::unique_ptr<ZdcNumberingScheme> ZdcSD::numberingScheme
private

Definition at line 44 of file ZdcSD.h.

Referenced by setDetUnitId(), and ZdcSD().

◆ showerLibrary

std::unique_ptr<ZdcShowerLibrary> ZdcSD::showerLibrary
private

Definition at line 43 of file ZdcSD.h.

Referenced by initRun(), and ZdcSD().

◆ thFibDir

double ZdcSD::thFibDir
private

Definition at line 41 of file ZdcSD.h.

◆ useShowerHits

bool ZdcSD::useShowerHits
private

Definition at line 39 of file ZdcSD.h.

Referenced by ProcessHits(), and ZdcSD().

◆ useShowerLibrary

bool ZdcSD::useShowerLibrary
private

Definition at line 39 of file ZdcSD.h.

Referenced by initRun(), and ZdcSD().

◆ verbosity

int ZdcSD::verbosity
private

Definition at line 38 of file ZdcSD.h.

Referenced by ZdcSD().

◆ zdcHitEnergyCut

double ZdcSD::zdcHitEnergyCut
private

Definition at line 42 of file ZdcSD.h.

Referenced by ZdcSD().