Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaGeneralProcess.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//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4GammaGeneralProcess
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 19.07.2018
37//
38// Modifications:
39//
40// Class Description:
41//
42// It is the gamma super process
43
44// -------------------------------------------------------------------
45//
46
47#ifndef G4GammaGeneralProcess_h
48#define G4GammaGeneralProcess_h 1
49
50#include "G4VEmProcess.hh"
51#include "globals.hh"
52#include "G4EmDataHandler.hh"
53
54class G4Step;
55class G4Track;
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
65{
66public:
67
68 explicit G4GammaGeneralProcess(const G4String& pname="GammaGeneralProc");
69
70 ~G4GammaGeneralProcess() override;
71
73
75
77
79
80 void ProcessDescription(std::ostream& outFile) const override;
81
82protected:
83
84 void InitialiseProcess(const G4ParticleDefinition*) override;
85
86public:
87
88 // Initialise for build of tables
89 void PreparePhysicsTable(const G4ParticleDefinition&) override;
90
91 // Build physics table during initialisation
92 void BuildPhysicsTable(const G4ParticleDefinition&) override;
93
94 // Called before tracking of each new G4Track
95 void StartTracking(G4Track*) override;
96
97 // implementation of virtual method, specific for G4GammaGeneralProcess
99 const G4Track& track,
100 G4double previousStepSize,
101 G4ForceCondition* condition) override;
102
103 // implementation of virtual method, specific for G4GammaGeneralProcess
104 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
105
106 // Store PhysicsTable in a file.
107 // Return false in case of failure at I/O
109 const G4String& directory,
110 G4bool ascii = false) override;
111
112 // Retrieve Physics from a file.
113 // (return true if the Physics Table can be build by using file)
114 // (return false if the process has no functionality or in case of failure)
115 // File name should is constructed as processName+particleName and the
116 // should be placed under the directory specifed by the argument.
118 const G4String& directory,
119 G4bool ascii) override;
120
121 // Return sub-process limiting current step
122 const G4VProcess* GetCreatorProcess() const override;
123 inline const G4VProcess* GetSelectedProcess() const;
124
125 // Temporary methods
126 const G4String& GetSubProcessName() const;
128
129 G4VEmProcess* GetEmProcess(const G4String& name) override;
130
131 inline G4HadronicProcess* GetGammaNuclear() const;
132
133 // hide copy constructor and assignment operator
136 (const G4GammaGeneralProcess &right) = delete;
137
138protected:
139
140 G4double GetMeanFreePath(const G4Track& track, G4double previousStepSize,
141 G4ForceCondition* condition) override;
142
143 inline G4double ComputeGeneralLambda(size_t idxe, size_t idxt);
144
145 inline G4double GetProbability(size_t idxt);
146
147 inline void SelectedProcess(const G4Step& step, G4VProcess* ptr);
148
149 inline void SelectEmProcess(const G4Step&, G4VEmProcess*);
150
151 void SelectHadProcess(const G4Track&, const G4Step&, G4HadronicProcess*);
152
153 // It returns the cross section per volume for energy/ material
155
156private:
157
158 G4bool RetrieveTable(G4VEmProcess*, const G4String& directory,
159 G4bool ascii);
160
161protected:
162
165
168
169
170private:
171 static G4EmDataHandler* theHandler;
172 static const size_t nTables = 15;
173 static G4bool theT[nTables];
174 static G4String nameT[nTables];
175
176 G4VEmProcess* thePhotoElectric = nullptr;
177 G4VEmProcess* theCompton = nullptr;
178 G4VEmProcess* theConversionEE = nullptr;
179 G4VEmProcess* theRayleigh = nullptr;
180 G4GammaConversionToMuons* theConversionMM = nullptr;
181
182 G4double minPEEnergy;
183 G4double minEEEnergy;
184 G4double minMMEnergy;
185 G4double peLambda = 0.0;
186
187 size_t nLowE = 40;
188 size_t nHighE = 50;
189 size_t idxEnergy = 0;
190};
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193
194inline G4double
196{
197 idxEnergy = idxe;
198 return factor*theHandler->GetVector(idxt, basedCoupleIndex)
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
212inline void
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
220
222{
224 SelectedProcess(step, proc);
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
228
230{
231 return selectedProc;
232}
233
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
236
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
243
244#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4PhysicsVector * GetVector(size_t itable, size_t ivec) const
void AddHadProcess(G4HadronicProcess *)
void SelectedProcess(const G4Step &step, G4VProcess *ptr)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void SelectHadProcess(const G4Track &, const G4Step &, G4HadronicProcess *)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double GetProbability(size_t idxt)
G4bool IsApplicable(const G4ParticleDefinition &) override
void AddEmProcess(G4VEmProcess *)
G4GammaGeneralProcess(G4GammaGeneralProcess &)=delete
void StartTracking(G4Track *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
const G4VProcess * GetSelectedProcess() const
void AddMMProcess(G4GammaConversionToMuons *)
G4HadronicProcess * GetGammaNuclear() const
G4double ComputeGeneralLambda(size_t idxe, size_t idxt)
void ProcessDescription(std::ostream &outFile) const override
const G4String & GetSubProcessName() const
G4GammaGeneralProcess(const G4String &pname="GammaGeneralProc")
G4VEmProcess * GetEmProcess(const G4String &name) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
const G4VProcess * GetCreatorProcess() const override
void InitialiseProcess(const G4ParticleDefinition *) override
void SelectEmProcess(const G4Step &, G4VEmProcess *)
G4HadronicProcess * theGammaNuclear
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const
void SetProcessDefinedStep(const G4VProcess *aValue)
G4StepPoint * GetPostStepPoint() const
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
size_t basedCoupleIndex
const G4MaterialCutsCouple * currentCouple
G4double preStepKinEnergy