Geant4  10.03.p01
 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 // $Id: G4HadronicProcess.hh 97172 2016-05-27 12:37:57Z gcosmo $
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"
53 #include "G4EnergyRangeManager.hh"
54 #include "G4Nucleus.hh"
55 #include "G4ReactionProduct.hh"
56 #include <vector>
59 
61 #include "G4HadronicProcessType.hh"
62 
63 class G4Track;
64 class G4Step;
65 class G4Element;
66 class G4ParticleChange;
68 
70 {
71 public:
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
86  const G4Element * elm,
87  const G4Material* mat = nullptr);
88 
89  // obsolete method to get cross section per element
90  inline
92  const G4Element * elm,
93  const G4Material* mat = nullptr)
94  { return GetElementCrossSection(part, elm, mat); }
95 
96  // generic PostStepDoIt recommended for all derived classes
97  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
98  const G4Step& aStep);
99 
100  // initialisation of physics tables and G4HadronicProcessStore
101  virtual void PreparePhysicsTable(const G4ParticleDefinition&);
102 
103  // build physics tables and print out the configuration of the process
104  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
105 
106  // dump physics tables
108  { theCrossSectionDataStore->DumpPhysicsTable(p); }
109 
110  // add cross section data set
111  inline void AddDataSet(G4VCrossSectionDataSet * aDataSet)
112  { theCrossSectionDataStore->AddDataSet(aDataSet);}
113 
114  // access to the list of hadronic interactions
115  std::vector<G4HadronicInteraction*>& GetHadronicInteractionList()
116  { return theEnergyRangeManager.GetHadronicInteractionList(); }
117 
118  // get inverse cross section per volume
119  G4double GetMeanFreePath(const G4Track &aTrack, G4double,
120  G4ForceCondition *);
121 
122  // access to the target nucleus
123  inline const G4Nucleus* GetTargetNucleus() const
124  { return &targetNucleus; }
125 
126  // G4ParticleDefinition* GetTargetDefinition();
127  inline const G4Isotope* GetTargetIsotope()
128  { return targetNucleus.GetIsotope(); }
129 
130  virtual void ProcessDescription(std::ostream& outFile) const;
131 
132 protected:
133 
134  // generic method to choose secondary generator
135  // recommended for all derived classes
137  const G4HadProjectile & aHadProjectile, G4Nucleus & aTargetNucleus,
138  G4Material* aMaterial, G4Element* anElement)
139  { return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile,
140  aTargetNucleus,
141  aMaterial,anElement);
142  }
143 
144  // access to the target nucleus
146  { return &targetNucleus; }
147 
148 public:
149 
150  void BiasCrossSectionByFactor(G4double aScale);
151 
152  // Integral option
153  inline void SetIntegral(G4bool val)
154  { useIntegralXS = val; }
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 
173  inline void MultiplyCrossSectionBy(G4double factor)
174  { aScaleFactor = factor; }
175 
176 protected:
177 
178  void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
179 
180  // access to the chosen generator
182  { return theInteraction; }
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,
195 
196  // Check 4-momentum balance
198 
199 private:
200  G4double XBiasSurvivalProbability();
201  G4double XBiasSecondaryWeight();
202 
203  // hide assignment operator as private
204  G4HadronicProcess& operator=(const G4HadronicProcess& right);
206 
207  // Set E/p conservation check levels from environment variables
208  void GetEnergyMomentumCheckEnvvars();
209 
210 protected:
211 
213 
215 
217 
218 private:
219 
220  G4EnergyRangeManager theEnergyRangeManager;
221 
222  G4HadronicInteraction* theInteraction;
223 
224  G4CrossSectionDataStore* theCrossSectionDataStore;
225 
226  G4HadronicProcessStore* theProcessStore;
227 
228  G4Nucleus targetNucleus;
229 
230  bool G4HadronicProcess_debug_flag;
231 
232  bool useIntegralXS;
233 
234  G4int nMatWarn;
235 
236  // Energy-momentum checking
237  std::pair<G4double, G4double> epCheckLevels;
238  G4bool levelsSetByProcess;
239 
240  std::vector<G4VLeadingParticleBiasing *> theBias;
241 
242  G4double theInitialNumberOfInteractionLength;
243 
244  G4double aScaleFactor;
245  G4bool xBiasOn;
246  G4double theLastCrossSection;
247 };
248 
249 #endif
250 
G4double G4ParticleHPJENDLHEData::G4double result
void AddDataSet(G4VCrossSectionDataSet *)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void BiasCrossSectionByFactor(G4double aScale)
G4double GetMicroscopicCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
G4HadronicInteraction * GetHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement) const
const char * p
Definition: xmltok.h:285
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
void SetIntegral(G4bool val)
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
G4HadProjectile thePro
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
int G4int
Definition: G4Types.hh:78
const G4Isotope * GetTargetIsotope()
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
G4HadronicProcessType
void RegisterMe(G4HadronicInteraction *a)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4ParticleChange * theTotalResult
const G4Nucleus * GetTargetNucleus() const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
bool G4bool
Definition: G4Types.hh:79
G4CrossSectionDataStore * GetCrossSectionDataStore()
Definition: G4Step.hh:76
virtual void ProcessDescription(std::ostream &outFile) const
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
void SetEpReportLevel(G4int level)
void DumpPhysicsTable(const G4ParticleDefinition &)
G4HadronicInteraction * GetHadronicInteraction() const
G4Nucleus * GetTargetNucleusPointer()
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:119
void DumpPhysicsTable(const G4ParticleDefinition &p)
void MultiplyCrossSectionBy(G4double factor)
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, G4Material *aMaterial, G4Element *anElement)
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
double G4double
Definition: G4Types.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4ForceCondition
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4double GetLastCrossSection()
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
G4ProcessType