Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32// G4HadronicProcess
33//
34// This is the top level Hadronic Process class
35// The inelastic, elastic, capture, and fission processes
36// should derive from this class
37//
38// original by H.P.Wellisch
39// J.L. Chuma, TRIUMF, 10-Mar-1997
40// Last modified: 04-Apr-1997
41// 19-May-2008 V.Ivanchenko cleanup and added comments
42// 05-Jul-2010 V.Ivanchenko cleanup commented lines
43// 28-Jul-2012 M.Maire add function GetTargetDefinition()
44// 14-Sep-2012 Inherit from RestDiscrete, use subtype code (now in ctor) to
45// configure base-class
46// 28-Sep-2012 M. Kelsey -- Undo inheritance change, keep new ctor
47
48#ifndef G4HadronicProcess_h
49#define G4HadronicProcess_h 1
50
51#include "globals.hh"
52#include "G4VDiscreteProcess.hh"
54#include "G4Nucleus.hh"
55#include "G4ReactionProduct.hh"
56#include <vector>
59
62
63class G4Track;
64class G4Step;
65class G4Element;
67
68
70{
71public:
72 G4HadronicProcess(const G4String& processName="Hadronic",
73 G4ProcessType procType=fHadronic);
74
75 // Preferred signature for subclasses, specifying their subtype here
76 G4HadronicProcess(const G4String& processName,
77 G4HadronicProcessType subType);
78
79 virtual ~G4HadronicProcess();
80
81 // register generator of secondaries
83
84 // get cross section per element
85 inline
87 const G4Element * elm,
88 const G4Material* mat = 0)
89 {
90 G4double x = theCrossSectionDataStore->GetCrossSection(part, elm, mat);
91 if(x < 0.0) { x = 0.0; }
92 return x;
93 }
94
95 // obsolete method to get cross section per element
96 inline
98 const G4Element * elm,
99 const G4Material* mat = 0)
100 { return GetElementCrossSection(part, elm, mat); }
101
102 // generic PostStepDoIt recommended for all derived classes
103 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
104 const G4Step& aStep);
105
106 // initialisation of physics tables and G4HadronicProcessStore
107 virtual void PreparePhysicsTable(const G4ParticleDefinition&);
108
109 // build physics tables and print out the configuration of the process
110 virtual void BuildPhysicsTable(const G4ParticleDefinition&);
111
112 // dump physics tables
114 { theCrossSectionDataStore->DumpPhysicsTable(p); }
115
116 // add cross section data set
117 inline void AddDataSet(G4VCrossSectionDataSet * aDataSet)
118 { theCrossSectionDataStore->AddDataSet(aDataSet);}
119
120 // access to the manager
122 { return &theEnergyRangeManager; }
123
124 // get inverse cross section per volume
127
128 // access to the target nucleus
129 inline const G4Nucleus* GetTargetNucleus() const
130 { return &targetNucleus; }
131
132 // G4ParticleDefinition* GetTargetDefinition();
134 { return targetNucleus.GetIsotope(); }
135
136 virtual void ProcessDescription(std::ostream& outFile) const;
137
138protected:
139
140 // generic method to choose secondary generator
141 // recommended for all derived classes
143 G4double kineticEnergy, G4Material* aMaterial, G4Element* anElement)
144 { return theEnergyRangeManager.GetHadronicInteraction(kineticEnergy,
145 aMaterial,anElement);
146 }
147
148 // access to the target nucleus
150 { return &targetNucleus; }
151
152public:
153
155
156 // Energy-momentum non-conservation limits and reporting
157 inline void SetEpReportLevel(G4int level)
158 { epReportLevel = level; }
159
160 inline void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
161 { epCheckLevels.first = relativeLevel;
162 epCheckLevels.second = absoluteLevel;
163 levelsSetByProcess = true;
164 }
165
166 inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const
167 { return epCheckLevels; }
168
169 // access to the cross section data store
171 {return theCrossSectionDataStore;}
172
174 { aScaleFactor = factor; }
175
176protected:
177
178 void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
179
180 // obsolete method will be removed
182 { return theEnergyRangeManager; }
183
184 // obsolete method will be removed
185 inline void SetEnergyRangeManager( const G4EnergyRangeManager &value )
186 { theEnergyRangeManager = value; }
187
188 // access to the chosen generator
190 { return theInteraction; }
191
192 // access to the cross section data set
194 { return theLastCrossSection; }
195
196 // fill result
197 void FillResult(G4HadFinalState* aR, const G4Track& aT);
198
199 // Check the result for catastrophic energy non-conservation
201 const G4Nucleus& targetNucleus,
202 G4HadFinalState* result) const;
203
204 // Check 4-momentum balance
206
207private:
208 G4double XBiasSurvivalProbability();
209 G4double XBiasSecondaryWeight();
210
211 // hide assignment operator as private
212 G4HadronicProcess& operator=(const G4HadronicProcess& right);
214
215 // Set E/p conservation check levels from environment variables
216 void GetEnergyMomentumCheckEnvvars();
217
218protected:
219
221
223
225
226private:
227
228 G4EnergyRangeManager theEnergyRangeManager;
229
230 G4HadronicInteraction* theInteraction;
231
232 G4CrossSectionDataStore* theCrossSectionDataStore;
233
234 G4Nucleus targetNucleus;
235
236 bool G4HadronicProcess_debug_flag;
237
238 // Energy-momentum checking
239 std::pair<G4double, G4double> epCheckLevels;
240 G4bool levelsSetByProcess;
241
242 std::vector<G4VLeadingParticleBiasing *> theBias;
243
244 G4double theInitialNumberOfInteractionLength;
245
246 G4double aScaleFactor;
247 G4bool xBiasOn;
248 G4double theLastCrossSection;
249};
250
251#endif
252
G4ForceCondition
G4HadronicProcessType
G4ProcessType
@ fHadronic
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void AddDataSet(G4VCrossSectionDataSet *)
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4HadronicInteraction * GetHadronicInteraction(const G4double kineticEnergy, const G4Material *aMaterial, const G4Element *anElement) const
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4HadProjectile thePro
const G4Nucleus * GetTargetNucleus() const
void SetEnergyRangeManager(const G4EnergyRangeManager &value)
void BiasCrossSectionByFactor(G4double aScale)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4Nucleus * GetTargetNucleusPointer()
void SetEpReportLevel(G4int level)
G4HadronicInteraction * GetHadronicInteraction() const
G4ParticleChange * theTotalResult
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4EnergyRangeManager * GetManagerPointer()
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4double GetLastCrossSection()
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4HadronicInteraction * ChooseHadronicInteraction(G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement)
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
G4double GetMicroscopicCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void DumpPhysicsTable(const G4ParticleDefinition &p)
void MultiplyCrossSectionBy(G4double factor)
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result) const
virtual void ProcessDescription(std::ostream &outFile) const
void RegisterMe(G4HadronicInteraction *a)
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
const G4Isotope * GetTargetIsotope()
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
const G4EnergyRangeManager & GetEnergyRangeManager() const
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:119
Definition: G4Step.hh:78
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76