Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QuasiCerenkov.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 G4Cerenkov, extended to support
28// offloading optical photon generation. Offloading can be enabled either
29// via the G4OpticalParameters::Instance()->SetCerenkovOffloadPhotons(true)
30// method or the UI command:
31//
32// /process/optical/cerenkov/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 G4CerenkovQuasiTrackInfo. 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 G4QuasiCerenkov_h
49#define G4QuasiCerenkov_h 1
50
51#include "globals.hh"
52#include "G4DynamicParticle.hh"
53#include "G4ForceCondition.hh"
54#include "G4GPILSelection.hh"
56#include "G4VProcess.hh"
57
58#include <map>
59
60class G4Material;
62class G4PhysicsTable;
63class G4Step;
64class G4Track;
66
68{
69 public:
70 explicit G4QuasiCerenkov(const G4String& processName = "QuasiCerenkov",
73
74 explicit G4QuasiCerenkov(const G4QuasiCerenkov& right);
75
76 G4QuasiCerenkov& operator=(const G4QuasiCerenkov& right) = delete;
77
78 G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
79 // Returns true -> 'is applicable', for all charged particles
80 // except short-lived particles.
81
82 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) override;
83 // Build table at a right time
84
85 void PreparePhysicsTable(const G4ParticleDefinition& part) override;
86 void Initialise();
87
89 // Returns the discrete step limit and sets the 'StronglyForced'
90 // condition for the DoIt to be invoked at every step.
91
93 G4ForceCondition*) override;
94 // Returns the discrete step limit and sets the 'StronglyForced'
95 // condition for the DoIt to be invoked at every step.
96
98 const G4Step& aStep) override;
99 // This is the method implementing the Cerenkov process.
100
101 // no operation in AtRestDoIt and AlongStepDoIt
103 const G4Track&, G4double, G4double, G4double&, G4GPILSelection*) override
104 {
105 return -1.0;
106 };
107
109 const G4Track&, G4ForceCondition*) override
110 {
111 return -1.0;
112 };
113
114 // no operation in AtRestDoIt and AlongStepDoIt
115 virtual G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override
116 {
117 return nullptr;
118 };
119
121 const G4Step&) override
122 {
123 return nullptr;
124 };
125
126 void SetTrackSecondariesFirst(const G4bool state);
127 // If set, the primary particle tracking is interrupted and any
128 // produced Cerenkov photons are tracked next. When all have
129 // been tracked, the tracking of the primary resumes.
130
132 // Returns the boolean flag for tracking secondaries first.
133
134 void SetMaxBetaChangePerStep(const G4double d);
135 // Set the maximum allowed change in beta = v/c in % (perCent) per step.
136
138 // Returns the maximum allowed change in beta = v/c in % (perCent)
139
140 void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
141 // Set the maximum number of Cerenkov photons allowed to be generated during
142 // a tracking step. This is an average ONLY; the actual number will vary
143 // around this average. If invoked, the maximum photon stack will roughly be
144 // of the size set. If not called, the step is not limited by the number of
145 // photons generated.
146
148 // Returns the maximum number of Cerenkov photons allowed to be
149 // generated during a tracking step.
150
151 void SetStackPhotons(const G4bool);
152 // Call by the user to set the flag for stacking the Cerenkov photons
153
154 G4bool GetStackPhotons() const;
155 // Return the boolean for whether or not the Cerenkov photons are stacked
156
157 void SetOffloadPhotons(const G4bool);
158 // Call by the user to set the flag for offloading the Cerenkov photons
159
161 // Return the boolean for whether or not the Cerenkov photons are offloaded
162
163 G4int GetNumPhotons() const;
164 // Returns the current number of scint. photons (after PostStepDoIt)
165
167 // Returns the address of the physics table.
168
169 void DumpPhysicsTable() const;
170 // Prints the physics table.
171
172 G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta,
173 const G4Material* aMaterial,
174 G4MaterialPropertyVector* Rindex) const;
175
176 void DumpInfo() const override {ProcessDescription(G4cout);};
177 void ProcessDescription(std::ostream& out) const override;
178
180 // sets verbosity
181
182 protected:
184 std::map<std::size_t, std::size_t> fIndexMPT;
185
186 private:
187 G4double fMaxBetaChange;
188
189 G4int fMaxPhotons;
190 G4int fNumPhotons;
191
192 G4bool fStackingFlag;
193 G4bool fOffloadingFlag;
194 G4bool fTrackSecondariesFirst;
195
196 G4int secID = -1; // creator modelID
197
198};
199
201{
202 return fTrackSecondariesFirst;
203}
204
206{
207 return fMaxBetaChange;
208}
209
210inline G4int G4QuasiCerenkov::GetMaxNumPhotonsPerStep() const { return fMaxPhotons; }
211
212inline G4bool G4QuasiCerenkov::GetStackPhotons() const { return fStackingFlag; }
213
214inline G4bool G4QuasiCerenkov::GetOffloadPhotons() const { return fOffloadingFlag; }
215
216inline G4int G4QuasiCerenkov::GetNumPhotons() const { return fNumPhotons; }
217
222
223#endif /* G4QuasiCerenkov_h */
G4ForceCondition
G4GPILSelection
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
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
void ProcessDescription(std::ostream &out) const override
void PreparePhysicsTable(const G4ParticleDefinition &part) override
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
void DumpPhysicsTable() const
G4QuasiCerenkov(const G4QuasiCerenkov &right)
G4PhysicsTable * thePhysicsTable
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
G4bool GetOffloadPhotons() const
void SetOffloadPhotons(const G4bool)
std::map< std::size_t, std::size_t > fIndexMPT
G4QuasiCerenkov(const G4String &processName="QuasiCerenkov", G4ProcessType type=fElectromagnetic)
G4QuasiCerenkov & operator=(const G4QuasiCerenkov &right)=delete
G4double GetMaxBetaChangePerStep() const
void DumpInfo() const override
G4int GetNumPhotons() const
G4PhysicsTable * GetPhysicsTable() const
G4bool GetStackPhotons() const
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
void SetTrackSecondariesFirst(const G4bool state)
void SetVerboseLevel(G4int)
void SetMaxBetaChangePerStep(const G4double d)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4bool GetTrackSecondariesFirst() const
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
void SetStackPhotons(const G4bool)
G4int GetMaxNumPhotonsPerStep() const
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46