Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronGeneralProcess.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: G4NeutronGeneralProcess
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 08.08.2022
37//
38// Modifications:
39//
40// Class Description:
41//
42// It is the neutron super process
43
44// -------------------------------------------------------------------
45//
46
47#ifndef G4NeutronGeneralProcess_h
48#define G4NeutronGeneralProcess_h 1
49
50#include "G4HadronicProcess.hh"
51#include "globals.hh"
52#include "G4HadDataHandler.hh"
53#include <vector>
54
55class G4Step;
56class G4Track;
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
65{
66public:
67
68 explicit G4NeutronGeneralProcess(const G4String& pname="NeutronGeneralProc");
69
70 ~G4NeutronGeneralProcess() override;
71
73
74 void ProcessDescription(std::ostream& outFile) const override;
75
76 // Initialise for build of tables
77 void PreparePhysicsTable(const G4ParticleDefinition&) override;
78
79 // Build physics table during initialisation
80 void BuildPhysicsTable(const G4ParticleDefinition&) override;
81
82 // Store internal tables after initialisation
84 const G4String& directory, G4bool ascii) override;
85
86 // Called before tracking of each new G4Track
87 void StartTracking(G4Track*) override;
88
89 // implementation of virtual method, specific for G4NeutronGeneralProcess
91 const G4Track& track,
92 G4double previousStepSize,
93 G4ForceCondition* condition) override;
94
95 // implementation of virtual method, specific for G4NeutronGeneralProcess
96 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
97
98 const G4VProcess* GetCreatorProcess() const override;
99
100 // Temporary method
101 const G4String& GetSubProcessName() const;
102
103 // Temporary method
105
109
110 inline const G4VProcess* GetSelectedProcess() const;
111
112 inline void SetTimeLimit(G4double val);
113
114 inline void SetMinEnergyLimit(G4double val);
115
116 // hide copy constructor and assignment operator
119 (const G4NeutronGeneralProcess &right) = delete;
120
121protected:
122
123 G4double GetMeanFreePath(const G4Track& track, G4double previousStepSize,
124 G4ForceCondition* condition) override;
125
126 inline G4double ComputeGeneralLambda(size_t idxe, size_t idxt);
127
128 inline G4double GetProbability(size_t idxt);
129
130 inline void SelectedProcess(const G4Step& step, G4HadronicProcess* ptr,
132
133private:
134
135 // partial cross section
136 G4double ComputeCrossSection(G4VCrossSectionDataSet*, const G4Material*,
137 G4double kinEnergy, G4double loge);
138
139 G4VCrossSectionDataSet* InitialisationXS(G4HadronicProcess*);
140
141 // total cross section
142 inline void CurrentCrossSection(const G4Track&);
143
144 static G4HadDataHandler* theHandler;
145 static const size_t nTables = 5;
146 static G4String nameT[nTables];
147
148 G4HadronicProcess* fInelastic = nullptr;
149 G4HadronicProcess* fElastic = nullptr;
150 G4HadronicProcess* fCapture = nullptr;
151 G4HadronicProcess* fSelectedProc = nullptr;
152
153 G4VCrossSectionDataSet* fInelasticXS = nullptr;
154 G4VCrossSectionDataSet* fElasticXS = nullptr;
155 G4VCrossSectionDataSet* fCaptureXS = nullptr;
156
157 G4CrossSectionDataStore* fXSSInelastic = nullptr;
158 G4CrossSectionDataStore* fXSSElastic = nullptr;
159 G4CrossSectionDataStore* fXSSCapture = nullptr;
160 G4CrossSectionDataStore* fCurrentXSS = nullptr;
161
162 const G4ParticleDefinition* fNeutron;
163 const G4Material* fCurrMat = nullptr;
164
165 G4double fMinEnergy;
166 G4double fMiddleEnergy;
167 G4double fMaxEnergy;
168 G4double fTimeLimit;
169 G4double fXSFactorInel = 1.0;
170 G4double fXSFactorEl = 1.0;
171 G4double fCurrE = 0.0;
172 G4double fCurrLogE = 0.0;
173 G4double fLambda = 0.0;
174
175 // number of bins per decade
176 std::size_t nLowE = 100;
177 std::size_t nHighE = 10;
178
179 std::size_t idxEnergy = 0;
180 std::size_t matIndex = 0;
181
182 G4bool isMaster = true;
183 std::vector<G4double> fXsec;
184};
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
188inline G4double
189G4NeutronGeneralProcess::ComputeGeneralLambda(std::size_t idxe, std::size_t idxt)
190{
191 idxEnergy = idxe;
192 return theHandler->GetVector(idxt, matIndex)
193 ->LogVectorValue(fCurrE, fCurrLogE);
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197
199{
200 return theHandler->GetVector(idxt, matIndex)
201 ->LogVectorValue(fCurrE, fCurrLogE);
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205
206inline void
210
211{
212 fSelectedProc = ptr;
213 fCurrentXSS = xs;
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
218
220{
221 return fSelectedProc;
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
225
226inline void G4NeutronGeneralProcess::CurrentCrossSection(const G4Track& track)
227{
228 G4double energy = track.GetKineticEnergy();
229 const G4Material* mat = track.GetMaterial();
230 if(mat != fCurrMat || energy != fCurrE) {
231 fCurrMat = mat;
232 matIndex = mat->GetIndex();
233 fCurrE = energy;
234 fCurrLogE = track.GetDynamicParticle()->GetLogKineticEnergy();
235 fLambda = (energy <= fMiddleEnergy) ? ComputeGeneralLambda(0, 0)
236 : ComputeGeneralLambda(1, 3);
237 currentInteractionLength = 1.0/fLambda;
238 }
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242
244{
245 fTimeLimit = val;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249
251{
252 fMinEnergy = val;
253}
254
255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256
257#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
G4double GetLogKineticEnergy() const
const G4PhysicsVector * GetVector(std::size_t itable, std::size_t ivec) const
size_t GetIndex() const
Definition: G4Material.hh:255
void SetCaptureProcess(G4HadronicProcess *)
void StartTracking(G4Track *) override
G4NeutronGeneralProcess(G4NeutronGeneralProcess &)=delete
void ProcessDescription(std::ostream &outFile) const override
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void SelectedProcess(const G4Step &step, G4HadronicProcess *ptr, G4CrossSectionDataStore *)
const G4VProcess * GetSelectedProcess() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void SetElasticProcess(G4HadronicProcess *)
G4double GetProbability(size_t idxt)
G4bool StorePhysicsTable(const G4ParticleDefinition *part, const G4String &directory, G4bool ascii) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double ComputeGeneralLambda(size_t idxe, size_t idxt)
const G4VProcess * GetCreatorProcess() const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4String & GetSubProcessName() const
void SetMinEnergyLimit(G4double val)
void SetInelasticProcess(G4HadronicProcess *)
G4bool IsApplicable(const G4ParticleDefinition &) override
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const
void SetProcessDefinedStep(const G4VProcess *aValue)
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
G4double currentInteractionLength
Definition: G4VProcess.hh:339