LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
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 //
55 // Scintillation Light Class Definition
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 //
78 
79 #ifndef OpFastScintillation_h
80 #define OpFastScintillation_h 1
81 
83 // Includes
85 
86 #include "Geant4/globals.hh"
87 #include "Geant4/templates.hh"
88 #include "Geant4/G4ThreeVector.hh"
89 #include "Geant4/G4ParticleMomentum.hh"
90 #include "Geant4/G4Step.hh"
91 #include "Geant4/G4VRestDiscreteProcess.hh"
92 #include "Geant4/G4OpticalPhoton.hh"
93 #include "Geant4/G4DynamicParticle.hh"
94 #include "Geant4/G4Material.hh"
95 #include "Geant4/G4PhysicsTable.hh"
96 #include "Geant4/G4MaterialPropertiesTable.hh"
97 #include "Geant4/G4PhysicsOrderedFreeVector.hh"
98 #include "Geant4/G4EmSaturation.hh"
99 
100 #include "fhiclcpp/ParameterSet.h"
101 #include "TF1.h"
102 
103 // Class Description:
104 // RestDiscrete Process - Generation of Scintillation Photons.
105 // Class inherits publicly from G4VRestDiscreteProcess.
106 // Class Description - End:
107 
109 // Class Definition
111 
112 namespace larg4{
113 
114 class OpFastScintillation : public G4VRestDiscreteProcess
115 {
116 
117 private:
118 
120  // Operators
122 
123  // OpFastScintillation& operator=(const OpFastScintillation &right);
124 
125 public: // Without description
126 
128  // Constructors and Destructor
130 
131  OpFastScintillation(const G4String& processName = "Scintillation", G4ProcessType type = fElectromagnetic);
133 
135 
137  // Methods
139 
140 public: // With description
141 
142  // OpFastScintillation Process has both PostStepDoIt (for energy
143  // deposition of particles in flight) and AtRestDoIt (for energy
144  // given to the medium by particles at rest)
145 
146  virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
147  // Returns true -> 'is applicable', for any particle type except
148  // for an 'opticalphoton' and for short-lived particles
149 
150  G4double GetMeanFreePath(const G4Track& aTrack,
151  G4double ,
152  G4ForceCondition* );
153  // Returns infinity; i. e. the process does not limit the step,
154  // but sets the 'StronglyForced' condition for the DoIt to be
155  // invoked at every step.
156 
157 
158  G4double GetMeanLifeTime(const G4Track& aTrack,
159  G4ForceCondition* );
160  // Returns infinity; i. e. the process does not limit the time,
161  // but sets the 'StronglyForced' condition for the DoIt to be
162  // invoked at every step.
163 
164  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
165  const G4Step& aStep);
166  virtual G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
167  const G4Step& aStep);
168 
169  // These are the methods implementing the scintillation process.
170 
171  void SetTrackSecondariesFirst(const G4bool state);
172  // If set, the primary particle tracking is interrupted and any
173  // produced scintillation photons are tracked next. When all
174  // have been tracked, the tracking of the primary resumes.
175 
176  void SetFiniteRiseTime(const G4bool state);
177  // If set, the OpFastScintillation process expects the user to have
178  // set the constant material property FAST/SLOWSCINTILLATIONRISETIME.
179 
180  G4bool GetTrackSecondariesFirst() const;
181  // Returns the boolean flag for tracking secondaries first.
182 
183  G4bool GetFiniteRiseTime() const;
184  // Returns the boolean flag for a finite scintillation rise time.
185 
186  void SetScintillationYieldFactor(const G4double yieldfactor);
187  // Called to set the scintillation photon yield factor, needed when
188  // the yield is different for different types of particles. This
189  // scales the yield obtained from the G4MaterialPropertiesTable.
190 
191  G4double GetScintillationYieldFactor() const;
192  // Returns the photon yield factor.
193 
194  void SetScintillationExcitationRatio(const G4double excitationratio);
195  // Called to set the scintillation exciation ratio, needed when
196  // the scintillation level excitation is different for different
197  // types of particles. This overwrites the YieldRatio obtained
198  // from the G4MaterialPropertiesTable.
199 
200  G4double GetScintillationExcitationRatio() const;
201  // Returns the scintillation level excitation ratio.
202 
203  G4PhysicsTable* GetFastIntegralTable() const;
204  // Returns the address of the fast scintillation integral table.
205 
206  G4PhysicsTable* GetSlowIntegralTable() const;
207  // Returns the address of the slow scintillation integral table.
208 
209  void AddSaturation(G4EmSaturation* sat) { emSaturation = sat; }
210  // Adds Birks Saturation to the process.
211 
212  void RemoveSaturation() { emSaturation = NULL; }
213  // Removes the Birks Saturation from the process.
214 
215  G4EmSaturation* GetSaturation() const { return emSaturation; }
216  // Returns the Birks Saturation.
217 
218  void SetScintillationByParticleType(const G4bool );
219  // Called by the user to set the scintillation yield as a function
220  // of energy deposited by particle type
221 
223  { return scintillationByParticleType; }
224  // Return the boolean that determines the method of scintillation
225  // production
226 
227  void DumpPhysicsTable() const;
228  // Prints the fast and slow scintillation integral tables.
229 
230  std::vector<double> GetVUVTime(double, int) const;
231  std::vector<double> GetVisibleTimeOnlyCathode(double, int) const;
232  // Update configuration parameters.
233 
234  //void reconfigure(const fhicl::ParameterSet& pset);
235 
236 protected:
237 
238  void BuildThePhysicsTable();
239  // It builds either the fast or slow scintillation integral table;
240  // or both.
241 
242 
243  bool RecordPhotonsProduced(const G4Step& aStep, double N);
244  // Note the production of N photons in at point xyz.
245  // pass on to generate detector response, etc.
246 
247 
249  // Class Data Members
251 
252  G4PhysicsTable* theSlowIntegralTable;
253  G4PhysicsTable* theFastIntegralTable;
254 
257 
258  G4double YieldFactor;
259 
260  G4double ExcitationRatio;
261 
263 
264 private:
265 
266  G4double single_exp(G4double t, G4double tau2) const;
267  G4double bi_exp(G4double t, G4double tau1, G4double tau2) const;
268 
269  G4double scint_time(const G4Step& aStep,
270  G4double ScintillationTime,
271  G4double ScintillationRiseTime) const;
272  std::vector<double> propagation_time(G4ThreeVector x0, int OpChannel, int NPhotons, bool Reflected=false) const;
273 
274  // emission time distribution when there is a finite rise time
275  G4double sample_time(G4double tau1, G4double tau2) const;
276 
277  // Facility for TPB emission energies
278  double reemission_energy() const;
279  std::map<double,double> tpbemission;
280  CLHEP::RandGeneral *rgen0;
281 
282  void average_position(G4Step const& aStep, double *xzyPos) const;
283 
284  G4EmSaturation* emSaturation;
285  // functions and parameters for the propagation time parametrization
287  float const* ReflT0s;
288 
289  TF1 const* functions_vuv[8];
290  TF1 const* functions_vis[5];
291  double fd_break;
292  double fd_max;
295  //double fGlobalTimeOffset;
296 
297  void ProcessStep( const G4Step& step);
298 
299  bool bPropagate;
300 
301 };
302 
303  double finter_d(double*, double*);
304  double LandauPlusExpoFinal(double*, double*);
305 
307 // Inline methods
309 
310 inline
311 G4bool OpFastScintillation::IsApplicable(const G4ParticleDefinition& aParticleType)
312 {
313  if (aParticleType.GetParticleName() == "opticalphoton") return false;
314  if (aParticleType.IsShortLived()) return false;
315 
316  return true;
317 }
318 
319 inline
321 {
322  fTrackSecondariesFirst = state;
323 }
324 
325 inline
327 {
328  fFiniteRiseTime = state;
329 }
330 
331 inline
333 {
334  return fTrackSecondariesFirst;
335 }
336 
337 inline
339 {
340  return fFiniteRiseTime;
341 }
342 
343 inline
344 void OpFastScintillation::SetScintillationYieldFactor(const G4double yieldfactor)
345 {
346  YieldFactor = yieldfactor;
347 }
348 
349 inline
351 {
352  return YieldFactor;
353 }
354 
355 inline
356 void OpFastScintillation::SetScintillationExcitationRatio(const G4double excitationratio)
357 {
358  ExcitationRatio = excitationratio;
359 }
360 
361 inline
363 {
364  return ExcitationRatio;
365 }
366 
367 inline
369 {
370  return theSlowIntegralTable;
371 }
372 
373 inline
375 {
376  return theFastIntegralTable;
377 }
378 
379 inline
381 {
382  if (theFastIntegralTable) {
383  G4int PhysicsTableSize = theFastIntegralTable->entries();
384  G4PhysicsOrderedFreeVector *v;
385 
386  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
387  {
388  v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i];
389  v->DumpValues();
390  }
391  }
392 
393  if (theSlowIntegralTable) {
394  G4int PhysicsTableSize = theSlowIntegralTable->entries();
395  G4PhysicsOrderedFreeVector *v;
396 
397  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
398  {
399  v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i];
400  v->DumpValues();
401  }
402  }
403 }
404 
405 inline
406 G4double OpFastScintillation::single_exp(G4double t, G4double tau2) const
407 {
408  return std::exp(-1.0*t/tau2)/tau2;
409 }
410 
411 inline
412 G4double OpFastScintillation::bi_exp(G4double t, G4double tau1, G4double tau2) const
413 {
414  return std::exp(-1.0*t/tau2)*(1-std::exp(-1.0*t/tau1))/tau2/tau2*(tau1+tau2);
415 }
416 
417 } //namespace
418 
419 #endif /* OpFastScintillation_h */
G4PhysicsTable * GetSlowIntegralTable() const
G4EmSaturation * GetSaturation() const
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:112
void average_position(G4Step const &aStep, double *xzyPos) const
G4double sample_time(G4double tau1, G4double tau2) const
G4double GetScintillationExcitationRatio() const
G4bool GetTrackSecondariesFirst() const
std::vector< double > propagation_time(G4ThreeVector x0, int OpChannel, int NPhotons, bool Reflected=false) const
Geant4 interface.
void ProcessStep(const G4Step &step)
double finter_d(double *x, double *par)
G4double scint_time(const G4Step &aStep, G4double ScintillationTime, G4double ScintillationRiseTime) const
std::map< double, double > tpbemission
OpFastScintillation(const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
G4double bi_exp(G4double t, G4double tau1, G4double tau2) const
virtual G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
void SetScintillationYieldFactor(const G4double yieldfactor)
void SetScintillationByParticleType(const G4bool)
bool RecordPhotonsProduced(const G4Step &aStep, double N)
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
void SetTrackSecondariesFirst(const G4bool state)
G4double single_exp(G4double t, G4double tau2) const
G4PhysicsTable * GetFastIntegralTable() const
std::vector< double > GetVisibleTimeOnlyCathode(double, int) const
G4bool GetScintillationByParticleType() const
void SetScintillationExcitationRatio(const G4double excitationratio)
std::vector< double > GetVUVTime(double, int) const
bool bPropagate
Whether propagation of photons is enabled.
void SetFiniteRiseTime(const G4bool state)
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
void AddSaturation(G4EmSaturation *sat)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
G4double GetScintillationYieldFactor() const
double LandauPlusExpoFinal(double *x, double *par)