Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QuasiScintillation.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25// Class Description:
26//
27// This class is a modified clone of G4Scintillation, extended to support
28// offloading optical photon generation. Offloading can be enabled either
29// via the G4OpticalParameters::Instance()->SetScintOffloadPhotons(true)
30// method or the UI command:
31//
32// /process/optical/scintillation/setOffloadPhotons true
33//
34// When offloading is enabled, the process generates a single secondary track
35// of type G4QuasiOpticalPhoton, along with associated metadata encapsulated
36// in G4ScintillationQuasiTrackInfo. This auxiliary track information is used
37// to generate optical photons at a later stage—typically during offloading.
38//
39// The intended workflow leverages G4VTrackingManager, which delegates these
40// secondary tracks to a dedicated G4ProcessManager for G4QuasiOpticalPhoton.
41// These tracks are then handled by a user-defined custom tracking manager,
42// independent of the default process managers used for other particles.
43//
44// The primary purpose of this class is to facilitate the transfer of essential
45// data for offloaded optical photon generation in heterogeneous computing
46// models
47
48#ifndef G4QuasiScintillation_h
49#define G4QuasiScintillation_h 1
50
51#include "globals.hh"
52#include "G4EmSaturation.hh"
53#include "G4OpticalPhoton.hh"
55
56#include <map>
57
58class G4PhysicsTable;
59class G4Step;
60class G4Track;
61
63{
64 public:
65 explicit G4QuasiScintillation(const G4String& procName = "QausiScintillation",
68
71
72 // G4QuasiScintillation Process has both PostStepDoIt (for energy
73 // deposition of particles in flight) and AtRestDoIt (for energy
74 // given to the medium by particles at rest)
75
76 G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
77 // Returns true -> 'is applicable', for any particle type except
78 // for an 'opticalphoton' and for short-lived particles
79
80 void ProcessDescription(std::ostream&) const override;
81 void DumpInfo() const override {ProcessDescription(G4cout);};
82
83 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) override;
84 // Build table at the right time
85
86 void PreparePhysicsTable(const G4ParticleDefinition& part) override;
87 void Initialise();
88
90 G4ForceCondition*) override;
91 // Returns infinity; i. e. the process does not limit the step,
92 // but sets the 'StronglyForced' condition for the DoIt to be
93 // invoked at every step.
94
95 G4double GetMeanLifeTime(const G4Track& aTrack, G4ForceCondition*) override;
96 // Returns infinity; i. e. the process does not limit the time,
97 // but sets the 'StronglyForced' condition for the DoIt to be
98 // invoked at every step.
99
101 const G4Step& aStep) override;
102 G4VParticleChange* AtRestDoIt(const G4Track& aTrack,
103 const G4Step& aStep) override;
104
106 const G4Step& aStep,
107 G4double& yield1,
108 G4double& yield2,
109 G4double& yield3,
110 G4double& timeconstant1,
111 G4double& timeconstant2,
112 G4double& timeconstant3);
113 // allow multiple time constants with scint by particle type
114 // Returns the number of scintillation photons calculated when
115 // scintillation depends on the particle type and energy
116 // deposited (includes nonlinear dependendency) and updates the
117 // yields for each channel
118
119 void SetTrackSecondariesFirst(const G4bool state);
120 // If set, the primary particle tracking is interrupted and any
121 // produced scintillation photons are tracked next. When all
122 // have been tracked, the tracking of the primary resumes.
123
125 // Returns the boolean flag for tracking secondaries first.
126
127 void SetFiniteRiseTime(const G4bool state);
128 // If set, the G4QuasiScintillation process expects the user to have
129 // set the constant material property SCINTILLATIONRISETIME{1,2,3}.
130
132 // Returns the boolean flag for a finite scintillation rise time.
133
135 // Returns the address of scintillation integral table #1.
136
138 // Returns the address of scintillation integral table #2.
139
141 // Returns the address of scintillation integral table #3.
142
143 void AddSaturation(G4EmSaturation* sat);
144 // Adds Birks Saturation to the process.
145
146 void RemoveSaturation();
147 // Removes the Birks Saturation from the process.
148
150 // Returns the Birks Saturation.
151
153 // Called by the user to set the scintillation yield as a function
154 // of energy deposited by particle type
155
157 // Return the boolean that determines the method of scintillation
158 // production
159
160 void SetScintillationTrackInfo(const G4bool trackType);
161 // Call by the user to set the G4ScintillationTrackInformation
162 // to scintillation photon track
163
165 // Return the boolean for whether or not the
166 // G4QuasiScintillationTrackInformation is set to the scint. photon track
167
168 void SetStackPhotons(const G4bool);
169 // Call by the user to set the flag for stacking the scint. photons
170
171 G4bool GetStackPhotons() const;
172 // Return the boolean for whether or not the scint. photons are stacked
173
174 void SetOffloadPhotons(const G4bool);
175 // Call by the user to set the flag for offloading the scint. photons
176
178 // Return the boolean for whether or not the scint. photons are offloaded
179
180 G4int GetNumPhotons() const;
181 // Returns the current number of scint. photons (after PostStepDoIt)
182
183 void DumpPhysicsTable() const;
184 // Prints the fast and slow scintillation integral tables.
185
187 // sets verbosity
188
189 private:
190 void BuildInverseCdfTable(const G4MaterialPropertyVector* MPV,
191 G4PhysicsFreeVector* vec) const;
192 // Build the inverse cumulative distribution function (C.D.F.) table
193 // for the scintillation photon energy spectrum
194
195 private:
196
197 G4PhysicsTable* fIntegralTable1;
198 G4PhysicsTable* fIntegralTable2;
199 G4PhysicsTable* fIntegralTable3;
200 std::map<std::size_t, std::size_t> fIndexMPT;
201
202 G4EmSaturation* fEmSaturation;
203 const G4ParticleDefinition* opticalphoton =
205
206 G4int fNumPhotons;
207
208 G4bool fScintillationByParticleType;
209 G4bool fScintillationTrackInfo;
210 G4bool fStackingFlag;
211 G4bool fOffloadingFlag;
212 G4bool fTrackSecondariesFirst;
213 G4bool fFiniteRiseTime;
214
215#ifdef G4DEBUG_SCINTILLATION
216 G4double ScintTrackEDep, ScintTrackYield;
217#endif
218 // emission time distribution when there is a finite rise time
219 G4double sample_time(G4double tau1, G4double tau2);
220
221 G4int secID = -1; // creator modelID
222 G4int fNumEnergyWarnings = 0;
223
224};
225
226////////////////////
227// Inline methods
228////////////////////
229
231{
232 return fTrackSecondariesFirst;
233}
234
236{
237 return fFiniteRiseTime;
238}
239
241{
242 return fIntegralTable1;
243}
244
246{
247 return fIntegralTable2;
248}
249
251{
252 return fIntegralTable3;
253}
254
256{
257 fEmSaturation = sat;
258}
259
261{
262 fEmSaturation = nullptr;
263}
264
266{
267 return fEmSaturation;
268}
269
271{
272 return fScintillationByParticleType;
273}
274
276{
277 return fScintillationTrackInfo;
278}
279
281{
282 return fStackingFlag;
283}
284
286{
287 return fOffloadingFlag;
288}
289
291{
292 return fNumPhotons;
293}
294
295#endif /* G4QuasiScintillation_h */
G4ForceCondition
G4PhysicsFreeVector G4MaterialPropertyVector
G4ProcessType
@ fElectromagnetic
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
static G4OpticalPhoton * OpticalPhotonDefinition()
G4QuasiScintillation(const G4String &procName="QausiScintillation", G4ProcessType type=fElectromagnetic)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4PhysicsTable * GetIntegralTable3() const
G4bool GetScintillationByParticleType() const
void DumpInfo() const override
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4bool GetTrackSecondariesFirst() const
void SetStackPhotons(const G4bool)
void SetScintillationTrackInfo(const G4bool trackType)
void ProcessDescription(std::ostream &) const override
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4QuasiScintillation(const G4QuasiScintillation &right)=delete
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
G4PhysicsTable * GetIntegralTable2() const
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *) override
void SetTrackSecondariesFirst(const G4bool state)
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep, G4double &yield1, G4double &yield2, G4double &yield3, G4double &timeconstant1, G4double &timeconstant2, G4double &timeconstant3)
G4EmSaturation * GetSaturation() const
G4PhysicsTable * GetIntegralTable1() const
void PreparePhysicsTable(const G4ParticleDefinition &part) override
void SetScintillationByParticleType(const G4bool)
void AddSaturation(G4EmSaturation *sat)
G4bool GetScintillationTrackInfo() const
void SetOffloadPhotons(const G4bool)
void SetFiniteRiseTime(const G4bool state)
G4QuasiScintillation & operator=(const G4QuasiScintillation &right)=delete
G4VRestDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)