Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4HadronicProcess.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// G4HadronicProcess
32//
33// This is the top level Hadronic Process class
34// The inelastic, elastic, capture, and fission processes
35// should derive from this class
36//
37// original by H.P.Wellisch
38// J.L. Chuma, TRIUMF, 10-Mar-1997
39// Last modified: 04-Apr-1997
40// 19-May-2008 V.Ivanchenko cleanup and added comments
41// 05-Jul-2010 V.Ivanchenko cleanup commented lines
42// 28-Jul-2012 M.Maire add function GetTargetDefinition()
43// 14-Sep-2012 Inherit from RestDiscrete, use subtype code (now in ctor) to
44// configure base-class
45// 28-Sep-2012 M. Kelsey -- Undo inheritance change, keep new ctor
46
47#ifndef G4HadronicProcess_h
48#define G4HadronicProcess_h 1
49
50#include "globals.hh"
51#include "G4VDiscreteProcess.hh"
53#include "G4Nucleus.hh"
54#include "G4ReactionProduct.hh"
57#include "G4Material.hh"
58#include "G4DynamicParticle.hh"
59#include "G4ThreeVector.hh"
60#include "G4HadXSTypes.hh"
61#include <vector>
62
63class G4Track;
64class G4Step;
65class G4Element;
72
74{
75public:
76 G4HadronicProcess(const G4String& processName="Hadronic",
77 G4ProcessType procType=fHadronic);
78
79 // Preferred signature for subclasses, specifying their subtype here
80 G4HadronicProcess(const G4String& processName,
81 G4HadronicProcessType subType);
82
83 ~G4HadronicProcess() override;
84
85 // register generator of secondaries
87
88 // get cross section per element
90 const G4Element * elm,
91 const G4Material* mat = nullptr);
92
93 // obsolete method to get cross section per element
94 inline
96 const G4Element * elm,
97 const G4Material* mat = nullptr);
98
99 // initialisation for a new track
100 void StartTracking(G4Track* track) override;
101
102 // compute step limit
104 G4double, G4ForceCondition*) override;
105
106 // generic PostStepDoIt recommended for all derived classes
107 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
108 const G4Step& aStep) override;
109
110 // initialisation of physics tables and G4HadronicProcessStore
111 void PreparePhysicsTable(const G4ParticleDefinition&) override;
112
113 // build physics tables and print out the configuration of the process
114 void BuildPhysicsTable(const G4ParticleDefinition&) override;
115
116 // dump physics tables
118
119 // add cross section data set
120 void AddDataSet(G4VCrossSectionDataSet * aDataSet);
121
122 // access to the list of hadronic interactions
123 std::vector<G4HadronicInteraction*>& GetHadronicInteractionList();
124
125 // access to an hadronic interaction by name
127
128 // access to the chosen generator
130
131 // get inverse cross section per volume
133 G4ForceCondition *) override;
134
135 // access to the target nucleus
136 inline const G4Nucleus* GetTargetNucleus() const;
137
139
140 inline const G4Isotope* GetTargetIsotope();
141
142 // methods needed for implementation of integral XS
144 const G4Material*,
145 const G4double kinEnergy);
146
147 inline G4HadXSType CrossSectionType() const;
148 inline void SetCrossSectionType(G4HadXSType val);
149
150 void ProcessDescription(std::ostream& outFile) const override;
151
152 // scale cross section
154 void MultiplyCrossSectionBy(G4double factor);
155 inline G4double CrossSectionFactor() const;
156
157 // Integral option
158 inline void SetIntegral(G4bool val);
159
160 // Energy-momentum non-conservation limits and reporting
161 inline void SetEpReportLevel(G4int level);
162 inline void SetEnergyMomentumCheckLevels(G4double relativeLevel,
163 G4double absoluteLevel);
164 inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const;
165
166 // access to the cross section data store
168
169 // access to the data for integral XS method
170 inline std::vector<G4TwoPeaksHadXS*>* TwoPeaksXS() const;
171 inline std::vector<G4double>* EnergyOfCrossSectionMax() const;
172
173 // hide assignment operator as private
176
177protected:
178
179 // generic method to choose secondary generator
180 // recommended for all derived classes
182 const G4HadProjectile & aHadProjectile, G4Nucleus& aTargetNucleus,
183 const G4Material* aMaterial, const G4Element* anElement);
184
185 // access to the cross section data set
187
188 // fill result
189 void FillResult(G4HadFinalState* aR, const G4Track& aT);
190
191 void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
192
193 // Check the result for catastrophic energy non-conservation
195 const G4Nucleus& targetNucleus,
196 G4HadFinalState* result);
197
198 // Check 4-momentum balance
200
201private:
202
203 void InitialiseLocal();
204 void UpdateCrossSectionAndMFP(const G4double kinEnergy);
205 void RecomputeXSandMFP(const G4double kinEnergy);
206
207 inline void DefineXSandMFP();
208 inline void ComputeXSandMFP();
209
210 G4double XBiasSurvivalProbability();
211 G4double XBiasSecondaryWeight();
212
213 // Set E/p conservation check levels from environment variables
214 void GetEnergyMomentumCheckEnvvars();
215
216protected:
217
219
222
228
230
231private:
232
233 G4EnergyRangeManager theEnergyRangeManager;
234 G4Nucleus targetNucleus;
235
236 G4HadronicInteraction* theInteraction = nullptr;
237 G4HadronicProcessStore* theProcessStore;
238 const G4HadronicProcess* masterProcess = nullptr;
239 const G4ParticleDefinition* firstParticle = nullptr;
240 const G4ParticleDefinition* currentParticle = nullptr;
241 const G4Material* currentMat = nullptr;
242 const G4DynamicParticle* fDynParticle = nullptr;
243
244 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
245 std::vector<G4TwoPeaksHadXS*>* fXSpeaks = nullptr;
246
247 G4double theMFP = DBL_MAX;
248 G4double minKinEnergy;
249
250 // counters
251 G4int nMatWarn = 0;
252 G4int nKaonWarn = 0;
253 G4int nICelectrons = 0;
254 G4int matIdx = 0;
255
256 // flags
257 G4bool levelsSetByProcess = false;
258 G4bool G4HadronicProcess_debug_flag = false;
259 G4bool useIntegralXS = true;
260 G4bool isMaster = true;
261
262 G4ThreeVector unitVector;
263
264 // Energy-momentum checking
265 std::pair<G4double, G4double> epCheckLevels;
266 std::vector<G4VLeadingParticleBiasing*> theBias;
267};
268
271 const G4Element * elm,
272 const G4Material* mat)
273{
274 return GetElementCrossSection(part, elm, mat);
275}
276
279{
280 return theInteraction;
281}
282
283inline const G4Nucleus*
285{
286 return &targetNucleus;
287}
288
290{
291 return targetNucleus.GetIsotope();
292}
293
294inline G4HadXSType
296{
297 return fXSType;
298}
299
300inline void
302{
303 fXSType = val;
304}
305
307{
308 return aScaleFactor;
309}
310
312{
313 useIntegralXS = val;
314}
315
317{
318 epReportLevel = level;
319}
320
321inline void
323 G4double absoluteLevel)
324{
325 epCheckLevels.first = relativeLevel;
326 epCheckLevels.second = absoluteLevel;
327 levelsSetByProcess = true;
328}
329
330inline std::pair<G4double, G4double>
332{
333 return epCheckLevels;
334}
335
338{
340}
341
342inline std::vector<G4TwoPeaksHadXS*>*
344{
345 return fXSpeaks;
346}
347
348inline std::vector<G4double>*
350{
351 return theEnergyOfCrossSectionMax;
352}
353
355ChooseHadronicInteraction(const G4HadProjectile& aHadProjectile,
356 G4Nucleus& aTargetNucleus,
357 const G4Material* aMaterial,
358 const G4Element* anElement)
359{
360 return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile,
361 aTargetNucleus,
362 aMaterial,anElement);
363}
364
366{
367 return &targetNucleus;
368}
369
371{
372 return theLastCrossSection;
373}
374
375inline void G4HadronicProcess::DefineXSandMFP()
376{
378 theCrossSectionDataStore->GetCrossSection(fDynParticle, currentMat);
379 theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
380}
381
382inline void G4HadronicProcess::ComputeXSandMFP()
383{
385 theCrossSectionDataStore->ComputeCrossSection(fDynParticle, currentMat);
386 theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
387}
388
389#endif
390
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
G4HadXSType
Definition: G4HadXSTypes.hh:44
@ fHadNoIntegral
Definition: G4HadXSTypes.hh:45
G4HadronicProcessType
G4ProcessType
@ fHadronic
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4HadronicInteraction * GetHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement) const
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4HadProjectile thePro
const G4Nucleus * GetTargetNucleus() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void ProcessDescription(std::ostream &outFile) const override
G4double ComputeCrossSection(const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
G4HadronicProcess & operator=(const G4HadronicProcess &right)=delete
void BiasCrossSectionByFactor(G4double aScale)
void StartTracking(G4Track *track) override
void SetCrossSectionType(G4HadXSType val)
G4Nucleus * GetTargetNucleusPointer()
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
void SetEpReportLevel(G4int level)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4HadronicInteraction * GetHadronicInteraction() const
G4ParticleChange * theTotalResult
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS() const
G4HadXSType CrossSectionType() const
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
~G4HadronicProcess() override
G4double GetLastCrossSection()
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * theCrossSectionDataStore
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4HadronicInteraction * GetHadronicModel(const G4String &)
G4double GetMicroscopicCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
G4HadronicProcess(const G4HadronicProcess &)=delete
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void DumpPhysicsTable(const G4ParticleDefinition &p)
G4double CrossSectionFactor() const
void MultiplyCrossSectionBy(G4double factor)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4ForceCondition *) override
void SetIntegral(G4bool val)
void RegisterMe(G4HadronicInteraction *a)
const G4Isotope * GetTargetIsotope()
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
std::vector< G4double > * EnergyOfCrossSectionMax() const
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:111
Definition: G4Step.hh:62
#define DBL_MAX
Definition: templates.hh:62