Geant4 10.7.0
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
51
52#include "G4VEmProcess.hh"
53#include "globals.hh"
54#include "G4EmDataHandler.hh"
55
56class G4Step;
57class G4Track;
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
67{
68public:
69
70 explicit G4GammaGeneralProcess();
71
72 virtual ~G4GammaGeneralProcess();
73
75
77
79
81
82 void ProcessDescription(std::ostream& outFile) const override;
83
84protected:
85
86 void InitialiseProcess(const G4ParticleDefinition*) override;
87
88public:
89
90 // Initialise for build of tables
91 void PreparePhysicsTable(const G4ParticleDefinition&) override;
92
93 // Build physics table during initialisation
94 void BuildPhysicsTable(const G4ParticleDefinition&) override;
95
96 // Called before tracking of each new G4Track
97 void StartTracking(G4Track*) override;
98
99 // implementation of virtual method, specific for G4GammaGeneralProcess
101 const G4Track& track,
102 G4double previousStepSize,
103 G4ForceCondition* condition) override;
104
105 // implementation of virtual method, specific for G4GammaGeneralProcess
106 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
107
108 // Store PhysicsTable in a file.
109 // Return false in case of failure at I/O
111 const G4String& directory,
112 G4bool ascii = false) override;
113
114 // Retrieve Physics from a file.
115 // (return true if the Physics Table can be build by using file)
116 // (return false if the process has no functionality or in case of failure)
117 // File name should is constructed as processName+particleName and the
118 // should be placed under the directory specifed by the argument.
120 const G4String& directory,
121 G4bool ascii) override;
122
123 const G4String& GetProcessName() const;
124
125 G4int GetProcessSubType() const;
126
127 G4VEmProcess* GetEmProcess(const G4String& name) override;
128
129 // hide copy constructor and assignment operator
132 (const G4GammaGeneralProcess &right) = delete;
133
134protected:
135
136 G4double GetMeanFreePath(const G4Track& track, G4double previousStepSize,
137 G4ForceCondition* condition) override;
138
139 inline G4double ComputeGeneralLambda(size_t idxe, size_t idxt);
140
141 inline G4double GetProbability(size_t idxt);
142
143 inline void SelectedProcess(const G4Step& track, const G4VProcess* ptr);
144
146 G4VEmProcess*);
147
150
151private:
152
153 // It returns the cross section per volume for energy/ material
154 G4double TotalCrossSectionPerVolume();
155
156protected:
157
160
161private:
162 static G4EmDataHandler* theHandler;
163 static const size_t nTables = 15;
164 static G4bool theT[nTables];
165 static G4String nameT[nTables];
166
167 G4VEmProcess* thePhotoElectric;
168 G4VEmProcess* theCompton;
169 G4VEmProcess* theConversionEE;
170 G4VEmProcess* theRayleigh;
171 G4GammaConversionToMuons* theConversionMM;
172
173 G4double minPEEnergy;
174 G4double minEEEnergy;
175 G4double minMMEnergy;
176 G4double peLambda;
177 G4double preStepLogE;
178 G4double factor;
179
180 size_t nLowE;
181 size_t nHighE;
182 size_t idxEnergy;
183 G4bool splineFlag;
184};
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
188inline G4double
190{
191 idxEnergy = idxe;
192 return factor*theHandler->GetVector(idxt, basedCoupleIndex)
193 ->LogVectorValue(preStepKinEnergy, preStepLogE);
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197
199{
200 return (theT[idxt]) ? theHandler->GetVector(idxt, basedCoupleIndex)
201 ->LogVectorValue(preStepKinEnergy, preStepLogE) : 1.0;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205
206inline void
208{
209 selectedProc = ptr;
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214
216 const G4Track& track, const G4Step& step, G4VEmProcess* proc)
217{
219 SelectedProcess(step, proc);
220 return proc->PostStepDoIt(track, step);
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224
225#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 *)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double GetProbability(size_t idxt)
G4bool IsApplicable(const G4ParticleDefinition &) override
G4VParticleChange * SampleHadSecondaries(const G4Track &, const G4Step &, G4HadronicProcess *)
void AddEmProcess(G4VEmProcess *)
G4GammaGeneralProcess(G4GammaGeneralProcess &)=delete
void StartTracking(G4Track *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
void AddMMProcess(G4GammaConversionToMuons *)
G4double ComputeGeneralLambda(size_t idxe, size_t idxt)
const G4String & GetProcessName() const
void ProcessDescription(std::ostream &outFile) const override
G4VEmProcess * GetEmProcess(const G4String &name) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
void SelectedProcess(const G4Step &track, const G4VProcess *ptr)
void InitialiseProcess(const G4ParticleDefinition *) override
const G4VProcess * selectedProc
G4HadronicProcess * theGammaNuclear
G4VParticleChange * SampleEmSecondaries(const G4Track &, const G4Step &, G4VEmProcess *)
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 theEnergy, const G4double theLogEnergy) const
void SetProcessDefinedStep(const G4VProcess *aValue)
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
size_t basedCoupleIndex
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
const G4MaterialCutsCouple * currentCouple
G4double preStepKinEnergy