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

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

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