source: trunk/source/processes/hadronic/management/include/G4HadronicProcess.hh@ 1331

Last change on this file since 1331 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 8.0 KB
Line 
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// This is the top level Hadronic Process class
29// The inelastic, elastic, capture, and fission processes
30// should derive from this class
31//
32// original by H.P.Wellisch
33// J.L. Chuma, TRIUMF, 10-Mar-1997
34// Last modified: 04-Apr-1997
35// 19-May-2008 V.Ivanchenko cleanup and added comments
36
37#ifndef G4HadronicProcess_h
38#define G4HadronicProcess_h 1
39
40#include "globals.hh"
41#include "G4VDiscreteProcess.hh"
42#include "G4EnergyRangeManager.hh"
43#include "G4Nucleus.hh"
44#include "G4ReactionProduct.hh"
45#include <vector>
46#include "G4VIsotopeProduction.hh"
47#include "G4IsoParticleChange.hh"
48#include "G4VCrossSectionDataSet.hh"
49#include "G4VLeadingParticleBiasing.hh"
50//#include "G4Delete.hh"
51#include "G4CrossSectionDataStore.hh"
52#include "G4HadronicProcessType.hh"
53
54class G4Track;
55class G4Step;
56class G4Element;
57class G4ParticleChange;
58
59
60class G4HadronicProcess : public G4VDiscreteProcess
61{
62public:
63
64 G4HadronicProcess(const G4String& processName = "Hadronic",
65 G4ProcessType aType = fHadronic);
66
67 virtual ~G4HadronicProcess();
68
69 // register generator of secondaries
70 void RegisterMe(G4HadronicInteraction* a);
71
72 // get cross section per element
73 virtual
74 G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle,
75 const G4Element *anElement,
76 G4double aTemp );
77
78 // generic PostStepDoIt recommended for all derived classes
79 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
80 const G4Step& aStep);
81
82 // initialisation of physics tables and G4HadronicProcessStore
83 virtual void PreparePhysicsTable(const G4ParticleDefinition&);
84
85 // build physics tables and print out the configuration of the process
86 virtual void BuildPhysicsTable(const G4ParticleDefinition&);
87
88 // dump physics tables
89 inline void DumpPhysicsTable(const G4ParticleDefinition& p)
90 { theCrossSectionDataStore->DumpPhysicsTable(p); }
91
92 // add cross section data set
93 inline void AddDataSet(G4VCrossSectionDataSet * aDataSet)
94 { theCrossSectionDataStore->AddDataSet(aDataSet);}
95
96 // access to the manager
97 inline G4EnergyRangeManager *GetManagerPointer()
98 { return &theEnergyRangeManager; }
99
100 // get inverse cross section per volume
101 G4double GetMeanFreePath(const G4Track &aTrack, G4double,
102 G4ForceCondition *);
103
104protected:
105
106 // reset number of interaction length and save
107 virtual void ResetNumberOfInteractionLengthLeft()
108 { G4VProcess::ResetNumberOfInteractionLengthLeft();
109 theInitialNumberOfInteractionLength =
110 G4VProcess::theNumberOfInteractionLengthLeft;
111 }
112
113 // generic method to choose secondary generator
114 // recommended for all derived classes
115 inline G4HadronicInteraction *ChooseHadronicInteraction(
116 G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement )
117 { return theEnergyRangeManager.GetHadronicInteraction(kineticEnergy,
118 aMaterial,anElement);
119 }
120
121public:
122
123 // Methods for isotope production
124 static void EnableIsotopeProductionGlobally();
125 static void DisableIsotopeProductionGlobally();
126
127 void EnableIsotopeCounting() {isoIsOnAnyway = 1;}
128 void DisableIsotopeCounting() {isoIsOnAnyway = -1;}
129
130 void RegisterIsotopeProductionModel(G4VIsotopeProduction * aModel)
131 { theProductionModels.push_back(aModel); }
132
133 static G4IsoParticleChange * GetIsotopeProductionInfo();
134
135 void BiasCrossSectionByFactor(G4double aScale);
136
137 // Energy-momentum non-conservation limits and reporting
138 inline void SetEpReportLevel(G4int level)
139 { epReportLevel = level; }
140
141 inline void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
142 { epCheckLevels.first = relativeLevel;
143 epCheckLevels.second = absoluteLevel;
144 levelsSetByProcess = true;
145 }
146
147 inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const
148 { return epCheckLevels; }
149
150protected:
151
152 // obsolete method will be removed
153 inline const G4EnergyRangeManager &GetEnergyRangeManager() const
154 { return theEnergyRangeManager; }
155
156 // obsolete method will be removed
157 inline void SetEnergyRangeManager( const G4EnergyRangeManager &value )
158 { theEnergyRangeManager = value; }
159
160 // access to the chosen generator
161 inline G4HadronicInteraction *GetHadronicInteraction()
162 { return theInteraction; }
163
164 // access to the cross section data store
165 inline G4CrossSectionDataStore* GetCrossSectionDataStore()
166 { return theCrossSectionDataStore; }
167
168 // access to the cross section data set
169 inline G4double GetLastCrossSection()
170 { return theLastCrossSection; }
171
172private:
173
174 void DumpState(const G4Track&, const G4String&);
175
176 void FillTotalResult(G4HadFinalState * aR, const G4Track & aT);
177
178 G4HadFinalState * DoIsotopeCounting(G4HadFinalState * aResult,
179 const G4Track & aTrack,
180 const G4Nucleus & aNucleus);
181
182 G4IsoResult * ExtractResidualNucleus(const G4Track & aTrack,
183 const G4Nucleus & aNucleus,
184 G4HadFinalState * aResult);
185
186 inline G4double GetTotalNumberOfInteractionLengthTraversed()
187 { return theInitialNumberOfInteractionLength
188 -G4VProcess::theNumberOfInteractionLengthLeft;
189 }
190
191 // inline void SetCrossSectionDataStore(G4CrossSectionDataStore* aDataStore)
192 // { theCrossSectionDataStore = aDataStore; }
193
194 G4double XBiasSurvivalProbability();
195 G4double XBiasSecondaryWeight();
196
197 void CheckEnergyMomentumConservation(const G4Track&, const G4Nucleus&);
198
199private:
200
201 G4EnergyRangeManager theEnergyRangeManager;
202
203 G4HadronicInteraction *theInteraction;
204
205 G4CrossSectionDataStore* theCrossSectionDataStore;
206
207 G4Nucleus targetNucleus;
208
209 G4HadronicProcess *dispatch;
210
211 bool G4HadronicProcess_debug_flag;
212
213 // Energy-momentum checking
214 G4int epReportLevel;
215 std::pair<G4double, G4double> epCheckLevels;
216 G4bool levelsSetByProcess;
217
218 // swiches for isotope production
219 static G4bool isoIsEnabled; // true or false; local swich overrides
220 G4int isoIsOnAnyway; // true(1), false(-1) or default(0)
221
222 G4IsoParticleChange theIsoPC;
223 std::vector<G4VIsotopeProduction *> theProductionModels;
224
225 std::vector<G4VLeadingParticleBiasing *> theBias;
226
227 static G4IsoParticleChange* theIsoResult;
228 static G4IsoParticleChange* theOldIsoResult;
229
230 G4ParticleChange* theTotalResult;
231
232 G4double theInitialNumberOfInteractionLength;
233
234 G4double aScaleFactor;
235 G4bool xBiasOn;
236 G4double theLastCrossSection;
237
238 G4int ModelingState;
239};
240
241#endif
242
Note: See TracBrowser for help on using the repository browser.