Geant4  10.01
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 86448 2014-11-12 09:48:41Z 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;
67 
69 {
70 public:
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  virtual ~G4HadronicProcess();
79 
80  // register generator of secondaries
82 
83  // get cross section per element
84  inline
86  const G4Element * elm,
87  const G4Material* mat = 0)
88  {
90  if(x < 0.0) { x = 0.0; }
91  return x;
92  }
93 
94  // obsolete method to get cross section per element
95  inline
97  const G4Element * elm,
98  const G4Material* mat = 0)
99  { return GetElementCrossSection(part, elm, mat); }
100 
101  // generic PostStepDoIt recommended for all derived classes
102  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
103  const G4Step& aStep);
104 
105  // initialisation of physics tables and G4HadronicProcessStore
106  virtual void PreparePhysicsTable(const G4ParticleDefinition&);
107 
108  // build physics tables and print out the configuration of the process
109  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
110 
111  // dump physics tables
114 
115  // add cross section data set
116  inline void AddDataSet(G4VCrossSectionDataSet * aDataSet)
117  { theCrossSectionDataStore->AddDataSet(aDataSet);}
118 
119  // access to the list of hadronic interactions
120  std::vector<G4HadronicInteraction*>& GetHadronicInteractionList()
122 
123  // get inverse cross section per volume
124  G4double GetMeanFreePath(const G4Track &aTrack, G4double,
125  G4ForceCondition *);
126 
127  // access to the target nucleus
128  inline const G4Nucleus* GetTargetNucleus() const
129  { return &targetNucleus; }
130 
131  // G4ParticleDefinition* GetTargetDefinition();
132  inline const G4Isotope* GetTargetIsotope()
133  { return targetNucleus.GetIsotope(); }
134 
135  virtual void ProcessDescription(std::ostream& outFile) const;
136 
137 protected:
138 
139  // generic method to choose secondary generator
140  // recommended for all derived classes
142  const G4HadProjectile & aHadProjectile, G4Nucleus & aTargetNucleus,
143  G4Material* aMaterial, G4Element* anElement)
144  { return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile,
145  aTargetNucleus,
146  aMaterial,anElement);
147  }
148 
149  // access to the target nucleus
151  { return &targetNucleus; }
152 
153 public:
154 
155  void BiasCrossSectionByFactor(G4double aScale);
156 
157  // Energy-momentum non-conservation limits and reporting
158  inline void SetEpReportLevel(G4int level)
159  { epReportLevel = level; }
160 
161  inline void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
162  { epCheckLevels.first = relativeLevel;
163  epCheckLevels.second = absoluteLevel;
164  levelsSetByProcess = true;
165  }
166 
167  inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const
168  { return epCheckLevels; }
169 
170  // access to the cross section data store
172  {return theCrossSectionDataStore;}
173 
175  { aScaleFactor = factor; }
176 
177 protected:
178 
179  void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
180 
181  // access to the chosen generator
183  { return theInteraction; }
184 
185  // access to the cross section data set
187  { return theLastCrossSection; }
188 
189  // fill result
190  void FillResult(G4HadFinalState* aR, const G4Track& aT);
191 
192  // Check the result for catastrophic energy non-conservation
194  const G4Nucleus& targetNucleus,
195  G4HadFinalState* result);
196 
197  // Check 4-momentum balance
199 
200 private:
203 
204  // hide assignment operator as private
207 
208  // Set E/p conservation check levels from environment variables
210 
211 protected:
212 
214 
216 
218 
219 private:
220 
222 
224 
226 
228 
230 
231  // Energy-momentum checking
232  std::pair<G4double, G4double> epCheckLevels;
234 
235  std::vector<G4VLeadingParticleBiasing *> theBias;
236 
238 
242 };
243 
244 #endif
245 
void AddDataSet(G4VCrossSectionDataSet *)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void BiasCrossSectionByFactor(G4double aScale)
G4HadronicInteraction * GetHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement) const
G4HadronicInteraction * theInteraction
G4EnergyRangeManager theEnergyRangeManager
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
G4double GetMicroscopicCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
G4HadProjectile thePro
G4double a
Definition: TRTMaterials.hh:39
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
int G4int
Definition: G4Types.hh:78
const G4Isotope * GetTargetIsotope()
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
G4HadronicProcessType
std::pair< G4double, G4double > epCheckLevels
void GetEnergyMomentumCheckEnvvars()
void RegisterMe(G4HadronicInteraction *a)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4CrossSectionDataStore * theCrossSectionDataStore
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4ParticleChange * theTotalResult
G4HadronicProcess & operator=(const G4HadronicProcess &right)
const G4Nucleus * GetTargetNucleus() const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4double XBiasSurvivalProbability()
bool G4bool
Definition: G4Types.hh:79
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4double theInitialNumberOfInteractionLength
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
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 &)
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
static const G4double factor
G4HadronicInteraction * GetHadronicInteraction() const
G4double XBiasSecondaryWeight()
G4Nucleus * GetTargetNucleusPointer()
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:119
void DumpPhysicsTable(const G4ParticleDefinition &p)
void MultiplyCrossSectionBy(G4double factor)
std::vector< G4VLeadingParticleBiasing * > theBias
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