G4S1Light.hh
Go to the documentation of this file.
1 #ifndef G4S1Light_h
2 #define G4S1Light_h 1
3 
4 #include "globals.hh"
5 #include "templates.hh"
6 #include "Randomize.hh"
7 #include "G4Poisson.hh"
8 #include "G4ThreeVector.hh"
9 #include "G4ParticleMomentum.hh"
10 #include "G4Step.hh"
11 #include "G4VRestDiscreteProcess.hh"
12 #include "G4OpticalPhoton.hh"
13 #include "G4DynamicParticle.hh"
14 #include "G4Material.hh"
15 #include "G4MaterialPropertiesTable.hh"
16 #include "G4PhysicsOrderedFreeVector.hh"
17 #include "G4ThermalElectron.hh"
18 
19 #include <G4SystemOfUnits.hh>
20 #include <G4PhysicalConstants.hh>
21 
22 //#include "LUXSimManager.hh"
23 
24 #define AVO 6.022e23 //Avogadro's number (#/mol)
25 #define EMASS 9.109e-31*CLHEP::kg
26 #define MillerDriftSpeed true
27 
28 #define GASGAP 0.25*CLHEP::cm //S2 generation region
29 #define BORDER 0*CLHEP::cm //liquid-gas border z-coordinate
30 
31 #define QE_EFF 1 //a base or maximum quantum efficiency
32 #define phe_per_e 1 //S2 gain for quick studies
33 
34 // different field regions, for gamma-X studies
35 #define WIN 0*CLHEP::mm //top Cu block (also, quartz window)
36 #define TOP 0 //top grid wires
37 #define ANE 0 //anode mesh
38 #define SRF 0 //liquid-gas interface
39 #define GAT 0 //gate grid
40 #define CTH 0 //cathode grid
41 #define BOT 0 //bottom PMT grid
42 #define PMT 0 //bottom Cu block and PMTs
43 
44 // this entire file is adapted from G4Scintillation.hh from Geant4.9.4
45 class G4S1Light : public G4VRestDiscreteProcess //class definition
46 {
47  // Class inherits publicly from G4VRestDiscreteProcess
48 private:
49 public: // constructor and destructor
50 
51  G4S1Light(const G4String& processName = "S1",
52  G4ProcessType type = fElectromagnetic);
53  ~G4S1Light();
54 
55 public: // methods, with descriptions
56  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
57  // Returns true -> 'is applicable', for any particle type except for an
58  // 'opticalphoton' and for short-lived particles
59 
60  G4double GetMeanFreePath(const G4Track& aTrack,
61  G4double ,
62  G4ForceCondition* );
63  // Returns infinity; i. e. the process does not limit the step, but
64  // sets the 'StronglyForced' condition for the DoIt to be invoked at
65  // every step.
66 
67  G4double GetMeanLifeTime(const G4Track& aTrack,
68  G4ForceCondition* );
69  // Returns infinity; i. e. the process does not limit the time, but
70  // sets the 'StronglyForced' condition for the DoIt to be invoked at
71  // every step.
72 
73  // For in-flight particles losing energy (or those stopped)
74  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
75  const G4Step& aStep);
76  G4VParticleChange* AtRestDoIt ( const G4Track& aTrack,
77  const G4Step& aStep);
78 
79  // These are the methods implementing the scintillation process.
80 
81  void SetTrackSecondariesFirst(const G4bool state);
82  // If set, the primary particle tracking is interrupted and any
83  // produced scintillation quanta are tracked next. When all have been
84  // tracked, the tracking of the primary resumes.
85 
86  G4bool GetTrackSecondariesFirst() const;
87  // Returns the boolean flag for tracking secondaries first.
88 
89  void SetScintillationYieldFactor(const G4double yieldfactor);
90  // Called to set the scintillation quantum yield factor, useful for
91  // shutting off scintillation entirely, or for producing a universal
92  // re-scaling to for example represent detector effects. Internally is
93  // used for Lindhard yield factor for NR. Default should be user-set
94  // to be 1 (for ER) in your simulation -- see NEST readme
95 
96  G4double GetScintillationYieldFactor() const;
97  // Returns the quantum (photon/electron) yield factor. See above.
98 
99  void SetScintillationExcitationRatio(const G4double excitationratio);
100  // Called to set the scintillation exciton-to-ion ratio, useful for
101  // when NumExcitons/NumIons is different for different types of parent
102  // particles. This overwrites the ExcitationRatio obtained from the
103  // G4MaterialPropertiesTable (e.g., 0.06 for LXe).
104 
105  G4double GetScintillationExcitationRatio() const;
106  // Returns the ratio of the number of excitons to ions. Read above.
107 
108 protected:
109  G4bool fTrackSecondariesFirst; // see above
110  //bools for tracking some special particle cases
112  G4double fKr83m;
113  G4double YieldFactor; // turns scint. on/off
114  G4double ExcitationRatio; // N_ex/N_i, the dimensionless ratio of
115  //initial excitons to ions
116 private:
117  //LUXSimManager *luxManager;
118 };
119 
120 ////////////////////
121 // Inline methods
122 ////////////////////
123 inline
124 G4bool G4S1Light::IsApplicable(const G4ParticleDefinition& aParticleType)
125 {
126  if (aParticleType.GetParticleName() == "opticalphoton") return false;
127  if (aParticleType.IsShortLived()) return false;
128  if (aParticleType.GetParticleName() == "thermalelectron") return false;
129  //if(abs(aParticleType.GetPDGEncoding())==2112 || //neutron (no E-dep.)
130  if(abs(aParticleType.GetPDGEncoding())==12 || //neutrinos (ditto)
131  abs(aParticleType.GetPDGEncoding())==14 ||
132  abs(aParticleType.GetPDGEncoding())==16) return false;
133 
134  return true;
135 }
136 
137 inline
138 void G4S1Light::SetTrackSecondariesFirst(const G4bool state)
139 {
140  fTrackSecondariesFirst = state;
141 }
142 
143 inline
145 {
146  return fTrackSecondariesFirst;
147 }
148 
149 inline
150 void G4S1Light::SetScintillationYieldFactor(const G4double yieldfactor)
151 {
152  YieldFactor = yieldfactor;
153 }
154 
155 inline
157 {
158  return YieldFactor;
159 }
160 
161 inline
162 void G4S1Light::SetScintillationExcitationRatio(const G4double excitationratio)
163 {
164  ExcitationRatio = excitationratio;
165 }
166 
167 inline
169 {
170  return ExcitationRatio;
171 }
172 
173 #endif /* G4S1Light_h */
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4S1Light.cc:1094
G4bool fTrackSecondariesFirst
Definition: G4S1Light.hh:109
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4S1Light.hh:138
G4double GetScintillationYieldFactor() const
Definition: G4S1Light.hh:156
G4double GetScintillationExcitationRatio() const
Definition: G4S1Light.hh:168
~G4S1Light()
Definition: G4S1Light.cc:97
G4double YieldFactor
Definition: G4S1Light.hh:113
T abs(T value)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4S1Light.hh:124
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
Definition: G4S1Light.cc:1110
void SetScintillationYieldFactor(const G4double yieldfactor)
Definition: G4S1Light.hh:150
G4bool GetTrackSecondariesFirst() const
Definition: G4S1Light.hh:144
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4S1Light.cc:100
G4bool fExcitedNucleus
Definition: G4S1Light.hh:111
G4double ExcitationRatio
Definition: G4S1Light.hh:114
void SetScintillationExcitationRatio(const G4double excitationratio)
Definition: G4S1Light.hh:162
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4S1Light.cc:109
G4bool fVeryHighEnergy
Definition: G4S1Light.hh:111
G4S1Light(const G4String &processName="S1", G4ProcessType type=fElectromagnetic)
Definition: G4S1Light.cc:82
G4double fKr83m
Definition: G4S1Light.hh:112
G4bool fMultipleScattering
Definition: G4S1Light.hh:111
G4bool fAlpha
Definition: G4S1Light.hh:111