OpFastScintillation.hh
Go to the documentation of this file.
1 // Class adapted for LArSoft by Ben Jones, MIT 10/10/12
2 //
3 // This class is a physics process based on the standard Geant4
4 // scintillation process.
5 //
6 // It has been stripped down and adapted to form the backbone of
7 // the LArG4 fast optical simulation. Photons, instead of being
8 // produced and added to the geant4 particle stack, are logged
9 // and used to predict the visibility of this step to each PMT in
10 // the detector.
11 //
12 // The photonvisibilityservice looks up the visibility of the relevant
13 // xyz point, and if a photon is detected at a given PMT, one OnePhoton
14 // object is logged in the OpDetPhotonTable
15 //
16 // At the end of the event, the OpDetPhotonTable is read out
17 // by LArG4, and detected photons are stored in the event.
18 //
19 // This process can be used alongside the standard Cerenkov process,
20 // which does step geant4 opticalphotons. Both the fast scintillation
21 // table and the geant4 sensitive detectors are read out by LArG4 to
22 // produce a combined SimPhoton collection.
23 //
24 //
25 //
26 // ********************************************************************
27 // * License and Disclaimer *
28 // * *
29 // * The Geant4 software is copyright of the Copyright Holders of *
30 // * the Geant4 Collaboration. It is provided under the terms and *
31 // * conditions of the Geant4 Software License, included in the file *
32 // * LICENSE and available at http://cern.ch/geant4/license . These *
33 // * include a list of copyright holders. *
34 // * *
35 // * Neither the authors of this software system, nor their employing *
36 // * institutes,nor the agencies providing financial support for this *
37 // * work make any representation or warranty, express or implied, *
38 // * regarding this software system or assume any liability for its *
39 // * use. Please see the license in the file LICENSE and URL above *
40 // * for the full disclaimer and the limitation of liability. *
41 // * *
42 // * This code implementation is the result of the scientific and *
43 // * technical work of the GEANT4 collaboration. *
44 // * By using, copying, modifying or distributing the software (or *
45 // * any work based on the software) you agree to acknowledge its *
46 // * use in resulting scientific publications, and indicate your *
47 // * acceptance of all terms of the Geant4 Software license. *
48 // ********************************************************************
49 //
50 //
51 // GEANT4 tag $Name: not supported by cvs2svn $
52 //
53 //
54 ////////////////////////////////////////////////////////////////////////
55 // Scintillation Light Class Definition
56 ////////////////////////////////////////////////////////////////////////
57 //
58 // File: OpFastScintillation.hh
59 // Description: Discrete Process - Generation of Scintillation Photons
60 // Version: 1.0
61 // Created: 1998-11-07
62 // Author: Peter Gumplinger
63 // Updated: 2010-10-20 Allow the scintillation yield to be a function
64 // of energy deposited by particle type
65 // Thanks to Zach Hartwig (Department of Nuclear
66 // Science and Engineeering - MIT)
67 // 2005-07-28 add G4ProcessType to constructor
68 // 2002-11-21 change to user G4Poisson for small MeanNumPotons
69 // 2002-11-07 allow for fast and slow scintillation
70 // 2002-11-05 make use of constant material properties
71 // 2002-05-16 changed to inherit from VRestDiscreteProcess
72 // 2002-05-09 changed IsApplicable method
73 // 1999-10-29 add method and class descriptors
74 //
75 // mail: gum@triumf.ca
76 //
77 ////////////////////////////////////////////////////////////////////////
78 
79 #ifndef OpFastScintillation_h
80 #define OpFastScintillation_h 1
81 
82 /////////////
83 // Includes
84 /////////////
85 
87 #include "larcoreobj/SimpleTypesAndConstants/geo_vectors.h" // geo::Point_t
88 #include "larsim/PhotonPropagation/PhotonVisibilityTypes.h" // phot::MappedT0s_t
89 
90 #include "Geant4/G4ForceCondition.hh"
91 #include "Geant4/G4ParticleDefinition.hh"
92 #include "Geant4/G4PhysicsOrderedFreeVector.hh"
93 #include "Geant4/G4PhysicsTable.hh"
94 #include "Geant4/G4ProcessType.hh"
95 #include "Geant4/G4String.hh"
96 #include "Geant4/G4ThreeVector.hh"
97 #include "Geant4/G4Types.hh"
98 #include "Geant4/G4VRestDiscreteProcess.hh"
99 
100 #include "TF1.h"
101 #include "TVector3.h"
102 
103 #include <memory> // std::unique_ptr
104 
105 class G4EmSaturation;
106 class G4Step;
107 class G4Track;
108 class G4VParticleChange;
109 namespace CLHEP {
110  class RandGeneral;
111 }
112 namespace geo {
113  class GeometryCore;
114 }
115 namespace phot {
116  class PhotonVisibilityService;
117 }
118 
119 // Class Description:
120 // RestDiscrete Process - Generation of Scintillation Photons.
121 // Class inherits publicly from G4VRestDiscreteProcess.
122 // Class Description - End:
123 
124 /////////////////////
125 // Class Definition
126 /////////////////////
127 
128 namespace larg4 {
129 
130  class OpFastScintillation : public G4VRestDiscreteProcess {
131 
132  private:
133  //////////////
134  // Operators
135  //////////////
136 
137  // OpFastScintillation& operator=(const OpFastScintillation &right);
138 
139  public: // Without description
140  ////////////////////////////////
141  // Constructors and Destructor
142  ////////////////////////////////
143 
144  OpFastScintillation(const G4String& processName = "Scintillation",
145  G4ProcessType type = fElectromagnetic);
146 
148 
149  ////////////
150  // Methods
151  ////////////
152 
153  public: // With description
154  // OpFastScintillation Process has both PostStepDoIt (for energy
155  // deposition of particles in flight) and AtRestDoIt (for energy
156  // given to the medium by particles at rest)
157 
158  virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
159  // Returns true -> 'is applicable', for any particle type except
160  // for an 'opticalphoton' and for short-lived particles
161 
162  G4double GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition*);
163  // Returns infinity; i. e. the process does not limit the step,
164  // but sets the 'StronglyForced' condition for the DoIt to be
165  // invoked at every step.
166 
167  G4double GetMeanLifeTime(const G4Track& aTrack, G4ForceCondition*);
168  // Returns infinity; i. e. the process does not limit the time,
169  // but sets the 'StronglyForced' condition for the DoIt to be
170  // invoked at every step.
171 
172  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
173  virtual G4VParticleChange* AtRestDoIt(const G4Track& aTrack, const G4Step& aStep);
174 
175  // These are the methods implementing the scintillation process.
176 
177  void SetTrackSecondariesFirst(const G4bool state);
178  // If set, the primary particle tracking is interrupted and any
179  // produced scintillation photons are tracked next. When all
180  // have been tracked, the tracking of the primary resumes.
181 
182  void SetFiniteRiseTime(const G4bool state);
183  // If set, the OpFastScintillation process expects the user to have
184  // set the constant material property FAST/SLOWSCINTILLATIONRISETIME.
185 
186  G4bool GetTrackSecondariesFirst() const;
187  // Returns the boolean flag for tracking secondaries first.
188 
189  G4bool GetFiniteRiseTime() const;
190  // Returns the boolean flag for a finite scintillation rise time.
191 
192  void SetScintillationYieldFactor(const G4double yieldfactor);
193  // Called to set the scintillation photon yield factor, needed when
194  // the yield is different for different types of particles. This
195  // scales the yield obtained from the G4MaterialPropertiesTable.
196 
197  G4double GetScintillationYieldFactor() const;
198  // Returns the photon yield factor.
199 
200  void SetScintillationExcitationRatio(const G4double excitationratio);
201  // Called to set the scintillation exciation ratio, needed when
202  // the scintillation level excitation is different for different
203  // types of particles. This overwrites the YieldRatio obtained
204  // from the G4MaterialPropertiesTable.
205 
206  G4double GetScintillationExcitationRatio() const;
207  // Returns the scintillation level excitation ratio.
208 
209  G4PhysicsTable* GetFastIntegralTable() const;
210  // Returns the address of the fast scintillation integral table.
211 
212  G4PhysicsTable* GetSlowIntegralTable() const;
213  // Returns the address of the slow scintillation integral table.
214 
215  void
216  AddSaturation(G4EmSaturation* sat)
217  {
218  emSaturation = sat;
219  }
220  // Adds Birks Saturation to the process.
221 
222  void
224  {
225  emSaturation = NULL;
226  }
227  // Removes the Birks Saturation from the process.
228 
229  G4EmSaturation*
231  {
232  return emSaturation;
233  }
234  // Returns the Birks Saturation.
235 
236  void SetScintillationByParticleType(const G4bool);
237  // Called by the user to set the scintillation yield as a function
238  // of energy deposited by particle type
239 
240  G4bool
242  {
243  return scintillationByParticleType;
244  }
245  // Return the boolean that determines the method of scintillation
246  // production
247 
248  void DumpPhysicsTable() const;
249  // Prints the fast and slow scintillation integral tables.
250 
251  /*std::vector<double> GetVUVTime(double, int);
252  std::vector<double> GetVisibleTimeOnlyCathode(double, int);*/
253  // old timings -- to be deleted
254 
255  void getVUVTimes(std::vector<double>& arrivalTimes, const double distance_in_cm, const size_t angle_bin);
256  void generateParam(const size_t index, const size_t angle_bin);
257  // Functions for vuv component Landau + Exponential timing parameterisation, updated method
258 
259  void getVISTimes(std::vector<double>& arrivalTimes, const TVector3 &ScintPoint, const TVector3 &OpDetPoint);
260  // Visible component timing parameterisation
261 
262  void detectedDirectHits(std::map<size_t, int>& DetectedNum,
263  const double Num,
264  geo::Point_t const& ScintPoint) const;
265  void detectedReflecHits(std::map<size_t, int>& ReflDetectedNum,
266  const double Num,
267  geo::Point_t const& ScintPoint) const;
268 
269  protected:
270  void BuildThePhysicsTable();
271  // It builds either the fast or slow scintillation integral table;
272  // or both.
273 
274  bool RecordPhotonsProduced(const G4Step& aStep, double N);
275  // Note the production of N photons in at point xyz.
276  // pass on to generate detector response, etc.
277 
278  ///////////////////////
279  // Class Data Members
280  ///////////////////////
281 
282  std::unique_ptr<G4PhysicsTable> theSlowIntegralTable;
283  std::unique_ptr<G4PhysicsTable> theFastIntegralTable;
284 
287 
288  G4double YieldFactor;
289 
290  G4double ExcitationRatio;
291 
293 
294  private:
295 
297  double h; // height
298  double w; // width
300  int type;
301  };
302 
303  /// Returns whether the semi-analytic visibility parametrization is being used.
304  bool usesSemiAnalyticModel() const;
305 
306  int VUVHits(const double Nphotons_created,
307  geo::Point_t const& ScintPoint,
308  OpticalDetector const& opDet) const;
309  // Calculates semi-analytic model number of hits for vuv component
310 
311  int VISHits(geo::Point_t const& ScintPoint,
312  OpticalDetector const& opDet,
313  const double cathode_hits_rec,
314  const std::array<double, 3> hotspot) const;
315  // Calculates semi-analytic model number of hits for visible component
316 
317  G4double single_exp(const G4double t, const G4double tau2) const;
318  G4double bi_exp(const G4double t, const G4double tau1, const G4double tau2) const;
319 
320  G4double scint_time(const G4Step& aStep,
321  G4double ScintillationTime,
322  G4double ScintillationRiseTime) const;
323  void propagationTime(std::vector<double>& arrival_time_dist,
324  G4ThreeVector x0,
325  const size_t OpChannel,
326  bool Reflected = false); //const;
327 
328  // emission time distribution when there is a finite rise time
329  G4double sample_time(const G4double tau1, const G4double tau2) const;
330 
331  // Facility for TPB emission energies
332  double reemission_energy() const;
333  std::map<double, double> tpbemission;
334  std::unique_ptr<CLHEP::RandGeneral> fTPBEm;
335 
336  void average_position(G4Step const& aStep, double* xzyPos) const;
337 
338  G4EmSaturation* emSaturation;
339  // functions and parameters for the propagation time parametrization
342 
343  /*TF1 const* functions_vuv[8];
344  TF1 const* functions_vis[5];
345  double fd_break;
346  double fd_max;
347  double ftf1_sampling_factor;
348  double ft0_max, ft0_break_point;*/
349 
350  //For new VUV time parametrization
351  double fstep_size, fmin_d, fmax_d, fvuv_vgroup_mean, fvuv_vgroup_max, finflexion_point_distance, fangle_bin_timing_vuv;
352  std::vector<std::vector<double>> fparameters[7];
353  // vector containing generated VUV timing parameterisations
354  std::vector<std::vector<TF1>> VUV_timing;
355  // vector containing min and max range VUV timing parameterisations are sampled to
356  std::vector<std::vector<double>> VUV_max;
357  std::vector<std::vector<double>> VUV_min;
358 
359  // For new VIS time parameterisation
360  double fvis_vmean, fangle_bin_timing_vis;
361  std::vector<double> fdistances_refl;
362  std::vector<double> fradial_distances_refl;
363  std::vector<std::vector<std::vector<double>>> fcut_off_pars;
364  std::vector<std::vector<std::vector<double>>> ftau_pars;
365 
366  struct Dims {
367  double h, w; // height, width
368  };
369 
370  // solid angle of rectangular aperture calculation functions
371  double Rectangle_SolidAngle(const double a, const double b, const double d) const;
372  double Rectangle_SolidAngle(Dims const& o, const std::array<double, 3> v) const;
373  // solid angle of circular aperture calculation functions
374  double Disk_SolidAngle(const double d, const double h, const double b) const;
375  // solid angle of a dome aperture calculation functions
376  double Omega_Dome_Model(const double distance, const double theta) const;
377 
378  // For VUV semi-analytic hits
379  // Gaisser-Hillas correction parameters for VUV Nhits estimation
380  G4double Gaisser_Hillas(const double x, const double* par) const;
382  // flat PDs
384  std::vector<std::vector<double>> fGHvuvpars_flat;
385  std::vector<double> fborder_corr_angulo_flat;
386  std::vector<std::vector<double>> fborder_corr_flat;
387  // dome PDs
389  std::vector<std::vector<double>> fGHvuvpars_dome;
390  std::vector<double> fborder_corr_angulo_dome;
391  std::vector<std::vector<double>> fborder_corr_dome;
392 
393  // For VIS semi-analytic hits
395  // correction parameters for VIS Nhits estimation
397  // flat PDs
398  std::vector<double> fvis_distances_x_flat;
399  std::vector<double> fvis_distances_r_flat;
400  std::vector<std::vector<std::vector<double>>> fvispars_flat;
401  // dome PDs
402  std::vector<double> fvis_distances_x_dome;
403  std::vector<double> fvis_distances_r_dome;
404  std::vector<std::vector<std::vector<double>>> fvispars_dome;
405 
406  // geometry properties
407  double fplane_depth, fcathode_zdimension, fcathode_ydimension;
408  TVector3 fcathode_centre;
409  std::vector<geo::BoxBoundedGeo> const fActiveVolumes;
410 
411 
412  // Optical detector properties for semi-analytic hits
413  // int foptical_detector_type; // unused
414  double fradius;
417  std::vector<geo::Point_t> fOpDetCenter;
418  std::vector<int> fOpDetType;
419  std::vector<double> fOpDetLength;
420  std::vector<double> fOpDetHeight;
421  //double fGlobalTimeOffset;
422 
423  void ProcessStep(const G4Step& step);
424 
425  bool const bPropagate; ///< Whether propagation of photons is enabled.
426 
427  /// Photon visibility service instance.
429 
430  /// Whether the semi-analytic model is being used for photon visibility.
431  bool const fUseNhitsModel = false;
432  /// Whether photon propagation is performed only from active volumes
433  bool const fOnlyActiveVolume = false;
434  /// Allows running even if light on cryostats `C:1` and higher is not supported.
435  /// Currently hard coded "no"
436  bool const fOnlyOneCryostat = false;
437  /// Whether the cathodes are fully opaque; currently hard coded "no".
438  bool const fOpaqueCathode = false;
439 
440  bool isOpDetInSameTPC(geo::Point_t const& ScintPoint, geo::Point_t const& OpDetPoint) const;
441  bool isScintInActiveVolume(geo::Point_t const& ScintPoint);
442  double interpolate(const std::vector<double>& xData,
443  const std::vector<double>& yData,
444  double x,
445  bool extrapolate,
446  size_t i = 0) const;
447  void interpolate3(std::array<double, 3>& inter,
448  const std::vector<double>& xData,
449  const std::vector<double>& yData1,
450  const std::vector<double>& yData2,
451  const std::vector<double>& yData3,
452  double x,
453  bool extrapolate) const;
454 
455  static std::vector<geo::BoxBoundedGeo> extractActiveVolumes(geo::GeometryCore const& geom);
456 
457  }; // class OpFastScintillation
458 
459  double finter_d(double*, double*);
460  double LandauPlusExpoFinal(double*, double*);
461  double model_close(double*, double*);
462  double model_far(double*, double*);
463 
464  static const size_t acos_bins = 2000000;
465  constexpr double acos_table(const double x);
466  double fast_acos(const double x);
467 
468  ////////////////////
469  // Inline methods
470  ////////////////////
471  inline G4bool
472  OpFastScintillation::IsApplicable(const G4ParticleDefinition& aParticleType)
473  {
474  if (aParticleType.GetParticleName() == "opticalphoton") return false;
475  if (aParticleType.IsShortLived()) return false;
476 
477  return true;
478  }
479 
480  inline void
481  OpFastScintillation::SetTrackSecondariesFirst(const G4bool state)
482  {
483  fTrackSecondariesFirst = state;
484  }
485 
486  inline void
487  OpFastScintillation::SetFiniteRiseTime(const G4bool state)
488  {
489  fFiniteRiseTime = state;
490  }
491 
492  inline G4bool
493  OpFastScintillation::GetTrackSecondariesFirst() const
494  {
495  return fTrackSecondariesFirst;
496  }
497 
498  inline G4bool
499  OpFastScintillation::GetFiniteRiseTime() const
500  {
501  return fFiniteRiseTime;
502  }
503 
504  inline void
505  OpFastScintillation::SetScintillationYieldFactor(const G4double yieldfactor)
506  {
507  YieldFactor = yieldfactor;
508  }
509 
510  inline G4double
511  OpFastScintillation::GetScintillationYieldFactor() const
512  {
513  return YieldFactor;
514  }
515 
516  inline void
517  OpFastScintillation::SetScintillationExcitationRatio(const G4double excitationratio)
518  {
519  ExcitationRatio = excitationratio;
520  }
521 
522  inline G4double
523  OpFastScintillation::GetScintillationExcitationRatio() const
524  {
525  return ExcitationRatio;
526  }
527 
528  inline G4PhysicsTable*
529  OpFastScintillation::GetSlowIntegralTable() const
530  {
531  return theSlowIntegralTable.get();
532  }
533 
534  inline G4PhysicsTable*
535  OpFastScintillation::GetFastIntegralTable() const
536  {
537  return theFastIntegralTable.get();
538  }
539 
540  inline void
541  OpFastScintillation::DumpPhysicsTable() const
542  {
543  if (theFastIntegralTable) {
544  G4int PhysicsTableSize = theFastIntegralTable->entries();
545  G4PhysicsOrderedFreeVector* v;
546  for (G4int i = 0; i < PhysicsTableSize; i++) {
547  v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i];
548  v->DumpValues();
549  }
550  }
551  if (theSlowIntegralTable) {
552  G4int PhysicsTableSize = theSlowIntegralTable->entries();
553  G4PhysicsOrderedFreeVector* v;
554  for (G4int i = 0; i < PhysicsTableSize; i++) {
555  v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i];
556  v->DumpValues();
557  }
558  }
559  }
560 
561  template <typename TReal>
562  inline constexpr double
563  dist(const TReal* x, const TReal* y, const unsigned int dimension)
564  {
565  double d = 0.;
566  for (unsigned int p = 0; p < dimension; ++p) {
567  d += (*(x + p) - *(y + p)) * (*(x + p) - *(y + p));
568  }
569  return std::sqrt(d);
570  }
571 
572  template <typename TVector3>
573  inline constexpr double
574  dist(const std::array<double, 3> x, const TVector3 y, const unsigned int dimension, const unsigned int start)
575  {
576  double d = 0.;
577  for (unsigned int p = start; p < dimension; ++p) {
578  d += (x[p] - y[p]) * (x[p] - y[p]);
579  }
580  return std::sqrt(d);
581  }
582 
583  // implements relative method - do not use for comparing with zero
584  // use this most of the time, tolerance needs to be meaningful in your context
585  template <typename TReal>
586  inline constexpr static bool
587  isApproximatelyEqual(TReal a, TReal b, TReal tolerance = std::numeric_limits<TReal>::epsilon())
588  {
589  TReal diff = std::fabs(a - b);
590  if (diff <= tolerance) return true;
591  if (diff < std::fmax(std::fabs(a), std::fabs(b)) * tolerance) return true;
592  return false;
593  }
594 
595  // supply tolerance that is meaningful in your context
596  // for example, default tolerance may not work if you are comparing double with float
597  template <typename TReal>
598  inline constexpr static bool
599  isApproximatelyZero(TReal a, TReal tolerance = std::numeric_limits<TReal>::epsilon())
600  {
601  if (std::fabs(a) <= tolerance) return true;
602  return false;
603  }
604 
605  // use this when you want to be on safe side
606  // for example, don't start rover unless signal is above 1
607  template <typename TReal>
608  inline constexpr static bool
609  isDefinitelyLessThan(TReal a, TReal b, TReal tolerance = std::numeric_limits<TReal>::epsilon())
610  {
611  TReal diff = a - b;
612  if (diff < tolerance) return true;
613  if (diff < std::fmax(std::fabs(a), std::fabs(b)) * tolerance) return true;
614  return false;
615  }
616 
617  template <typename TReal>
618  inline constexpr static bool
619  isDefinitelyGreaterThan(TReal a, TReal b, TReal tolerance = std::numeric_limits<TReal>::epsilon())
620  {
621  TReal diff = a - b;
622  if (diff > tolerance) return true;
623  if (diff > std::fmax(std::fabs(a), std::fabs(b)) * tolerance) return true;
624  return false;
625  }
626 
627  // template<typename Function, typename... Args>
628  // auto OpFastScintillation::invoke_memoized(Function function, Args... args)
629  // {
630  // using key_type = std::tuple<Args...>;
631  // using value_type = std::invoke_result_t<Function, Args...>;
632  // static_assert(! std::is_same_v<Function, std::function<value_type(Args...)>>,
633  // "cannot memoize on std::function (use a lambda instead)");
634  // static_assert(! std::is_same_v<Function, value_type(*)(Args...)>,
635  // "cannot memoize on function pointer (use a lambda instead)");
636  // static std::mutex mutex;
637  // static std::map<key_type, value_type> cache;
638  // auto key = std::tuple(args...);
639  // auto lock = std::lock_guard<std::mutex>(mutex);
640  // if (cache.count(key)){
641  // return cache[key];
642  // }
643  // return cache[key] = std::apply(function, key);
644  // }
645 
646 } // namespace larg4
647 
648 #endif /* OpFastScintillation_h */
Index OpChannel(Index detNum, Index channel)
std::vector< std::vector< std::vector< double > > > fvispars_dome
G4EmSaturation * GetSaturation() const
std::vector< geo::Point_t > fOpDetCenter
std::vector< double > fOpDetHeight
auto const tolerance
std::map< double, double > tpbemission
std::vector< std::vector< std::vector< double > > > fvispars_flat
Geant4 interface.
constexpr double dist(const std::array< double, 3 > x, const TVector3 y, const unsigned int dimension, const unsigned int start)
std::vector< double > fvis_distances_r_flat
G4double tau1[100]
Definition: G4S2Light.cc:61
double finter_d(double *x, double *par)
std::vector< double > fdistances_refl
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
double model_far(double *x, double *par)
std::vector< std::vector< TF1 > > VUV_timing
double model_close(double *x, double *par)
std::vector< std::vector< double > > fGHvuvpars_flat
std::vector< std::vector< double > > VUV_min
std::vector< std::vector< double > > fGHvuvpars_dome
std::vector< double > fborder_corr_angulo_dome
constexpr unsigned int dimension()
bool const bPropagate
Whether propagation of photons is enabled.
std::vector< std::vector< double > > fborder_corr_dome
std::vector< std::vector< double > > fborder_corr_flat
std::vector< double > fvis_distances_x_flat
constexpr double acos_table(const double x)
std::vector< double > fborder_corr_angulo_flat
std::vector< double > fOpDetLength
const double a
std::unique_ptr< CLHEP::RandGeneral > fTPBEm
std::vector< geo::BoxBoundedGeo > const fActiveVolumes
std::vector< std::vector< std::vector< double > > > ftau_pars
p
Definition: test.py:223
double distance(double x1, double y1, double z1, double x2, double y2, double z2)
G4bool GetScintillationByParticleType() const
Description of geometry of one entire detector.
static constexpr bool isDefinitelyGreaterThan(TReal a, TReal b, TReal tolerance=std::numeric_limits< TReal >::epsilon())
Provides a base class aware of world box coordinates.
std::vector< double > fradial_distances_refl
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
std::vector< double > fvis_distances_r_dome
General LArSoft Utilities.
std::vector< std::vector< std::vector< double > > > fcut_off_pars
A container for photon visibility mapping data.
static constexpr bool isApproximatelyEqual(TReal a, TReal b, TReal tolerance=std::numeric_limits< TReal >::epsilon())
std::vector< double > fvis_distances_x_dome
Declaration of types related to photon visibility.
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
std::vector< std::vector< double > > VUV_max
static bool * b
Definition: config.cpp:1043
Definitions of geometry vector data types.
list x
Definition: train.py:276
void AddSaturation(G4EmSaturation *sat)
static const size_t acos_bins
static constexpr bool isDefinitelyLessThan(TReal a, TReal b, TReal tolerance=std::numeric_limits< TReal >::epsilon())
phot::MappedFunctions_t ParPropTimeTF1
LArSoft geometry interface.
Definition: ChannelGeo.h:16
static constexpr bool isApproximatelyZero(TReal a, TReal tolerance=std::numeric_limits< TReal >::epsilon())
double LandauPlusExpoFinal(double *x, double *par)
double fast_acos(double x)
std::unique_ptr< G4PhysicsTable > theFastIntegralTable