Geant4-11
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 <vector>
58
59class G4Track;
60class G4Step;
61class G4Element;
67
69{
70public:
71 G4HadronicProcess(const G4String& processName="Hadronic",
72 G4ProcessType procType=fHadronic);
73
74 // Preferred signature for subclasses, specifying their subtype here
75 G4HadronicProcess(const G4String& processName,
76 G4HadronicProcessType subType);
77
78 ~G4HadronicProcess() override;
79
80 // register generator of secondaries
82
83 // get cross section per element
85 const G4Element * elm,
86 const G4Material* mat = nullptr);
87
88 // obsolete method to get cross section per element
89 inline
91 const G4Element * elm,
92 const G4Material* mat = nullptr)
93 { return GetElementCrossSection(part, elm, mat); }
94
95 // generic PostStepDoIt recommended for all derived classes
97 const G4Step& aStep) override;
98
99 // initialisation of physics tables and G4HadronicProcessStore
100 void PreparePhysicsTable(const G4ParticleDefinition&) override;
101
102 // build physics tables and print out the configuration of the process
103 void BuildPhysicsTable(const G4ParticleDefinition&) override;
104
105 // dump physics tables
107
108 // add cross section data set
109 void AddDataSet(G4VCrossSectionDataSet * aDataSet);
110
111 // access to the list of hadronic interactions
112 std::vector<G4HadronicInteraction*>& GetHadronicInteractionList();
113
114 // access to an hadronic interaction by name
116
117 // access to the chosen generator
119 { return theInteraction; }
120
121 // get inverse cross section per volume
123 G4ForceCondition *) override;
124
125 // access to the target nucleus
126 inline const G4Nucleus* GetTargetNucleus() const
127 { return &targetNucleus; }
128
129 // G4ParticleDefinition* GetTargetDefinition();
131 { return targetNucleus.GetIsotope(); }
132
133 void ProcessDescription(std::ostream& outFile) const override;
134
135protected:
136
137 // generic method to choose secondary generator
138 // recommended for all derived classes
140 const G4HadProjectile & aHadProjectile, G4Nucleus& aTargetNucleus,
141 const G4Material* aMaterial, const G4Element* anElement)
142 { return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile,
143 aTargetNucleus,
144 aMaterial,anElement);
145 }
146
147 // access to the target nucleus
149 { return &targetNucleus; }
150
151public:
152
153 // scale cross section
155 void MultiplyCrossSectionBy(G4double factor);
157 { return aScaleFactor; }
158
159 // Integral option
160 inline void SetIntegral(G4bool val)
161 { useIntegralXS = val; }
162
163 // Energy-momentum non-conservation limits and reporting
164 inline void SetEpReportLevel(G4int level)
165 { epReportLevel = level; }
166
167 inline void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
168 { epCheckLevels.first = relativeLevel;
169 epCheckLevels.second = absoluteLevel;
170 levelsSetByProcess = true;
171 }
172
173 inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const
174 { return epCheckLevels; }
175
176 // access to the cross section data store
179
180protected:
181
182 void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
183
184 // access to the cross section data set
186 { return theLastCrossSection; }
187
188 // fill result
189 void FillResult(G4HadFinalState* aR, const G4Track& aT);
190
191 // Check the result for catastrophic energy non-conservation
193 const G4Nucleus& targetNucleus,
194 G4HadFinalState* result);
195
196 // Check 4-momentum balance
198
199private:
200
201 void InitialiseLocal();
202
205
206 // hide assignment operator as private
209
210 // Set E/p conservation check levels from environment variables
212
213protected:
214
216
218
220
222
223private:
224
226
228
230
232
234
236
238
241
243
244 // Energy-momentum checking
245 std::pair<G4double, G4double> epCheckLevels;
247
248 std::vector<G4VLeadingParticleBiasing*> theBias;
249
251
254};
255
256#endif
257
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
G4HadronicProcessType
G4ProcessType
@ fHadronic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
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
void BiasCrossSectionByFactor(G4double aScale)
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)
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
~G4HadronicProcess() override
G4HadronicProcess & operator=(const G4HadronicProcess &right)
G4double GetLastCrossSection()
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * theCrossSectionDataStore
G4double theInitialNumberOfInteractionLength
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4HadronicInteraction * GetHadronicModel(const G4String &)
std::vector< G4VLeadingParticleBiasing * > theBias
G4EnergyRangeManager theEnergyRangeManager
G4double GetMicroscopicCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
G4double XBiasSecondaryWeight()
G4HadronicProcessStore * theProcessStore
G4HadronicInteraction * theInteraction
G4double XBiasSurvivalProbability()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void GetEnergyMomentumCheckEnvvars()
void DumpPhysicsTable(const G4ParticleDefinition &p)
G4double CrossSectionFactor() const
void MultiplyCrossSectionBy(G4double factor)
std::pair< G4double, G4double > epCheckLevels
void SetIntegral(G4bool val)
void RegisterMe(G4HadronicInteraction *a)
G4HadronicProcess(const G4HadronicProcess &)
const G4Isotope * GetTargetIsotope()
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:111
Definition: G4Step.hh:62