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

Last change on this file since 900 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 8.2 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 // This is an abstract base class, since the pure virtual function
32 // PostStepDoIt has not been defined yet.
33 // Note: there is no .cc file
34 //
35 // original by H.P.Wellisch
36 // J.L. Chuma, TRIUMF, 10-Mar-1997
37 // Last modified: 04-Apr-1997
38
39#ifndef G4HadronicProcess_h
40#define G4HadronicProcess_h 1
41
42#include "globals.hh"
43#include "G4VDiscreteProcess.hh"
44#include "G4EnergyRangeManager.hh"
45#include "G4Nucleus.hh"
46#include "G4ReactionProduct.hh"
47#include <vector>
48#include "G4VIsotopeProduction.hh"
49#include "G4IsoParticleChange.hh"
50#include "G4VCrossSectionDataSet.hh"
51#include "G4VLeadingParticleBiasing.hh"
52#include "G4Delete.hh"
53#include "G4CrossSectionDataStore.hh"
54#include "G4HadronicException.hh"
55
56
57class G4Track;
58class G4Step;
59class G4Element;
60class G4ParticleChange;
61
62 class G4HadronicProcess : public G4VDiscreteProcess
63 {
64 public:
65
66 G4HadronicProcess( const G4String &processName = "Hadronic",
67 G4ProcessType aType = fHadronic );
68 virtual ~G4HadronicProcess();
69
70 void RegisterMe( G4HadronicInteraction *a );
71
72 void AddDataSet(G4VCrossSectionDataSet * aDataSet)
73 {
74 theCrossSectionDataStore->AddDataSet(aDataSet);
75 }
76
77 virtual G4VParticleChange *PostStepDoIt( const G4Track &aTrack,
78 const G4Step &aStep ) = 0;
79
80 virtual
81 G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle,
82 const G4Element *anElement,
83 G4double aTemp ) = 0;
84
85 G4double GetMeanFreePath(const G4Track &aTrack, G4double,
86 G4ForceCondition *);
87
88 // Set methods for isotope production
89
90 static void EnableIsotopeProductionGlobally();
91 static void DisableIsotopeProductionGlobally();
92
93 void EnableIsotopeCounting() {isoIsOnAnyway = 1;}
94 void DisableIsotopeCounting() {isoIsOnAnyway = -1;}
95
96 void RegisterIsotopeProductionModel(G4VIsotopeProduction * aModel)
97 { theProductionModels.push_back(aModel); }
98
99 static G4IsoParticleChange * GetIsotopeProductionInfo()
100 {
101 G4IsoParticleChange * anIsoResult = theIsoResult;
102 if(theIsoResult) theOldIsoResult = theIsoResult;
103 theIsoResult = 0;
104 return anIsoResult;
105 }
106
107 void BiasCrossSectionByFactor(G4double aScale)
108 {
109 xBiasOn = true;
110 aScaleFactor = aScale;
111 G4String it = GetProcessName();
112 if( (it != "PhotonInelastic") &&
113 (it != "ElectroNuclear") &&
114 (it != "PositronNuclear") )
115 {
116 G4Exception("G4HadronicProcess", "007", FatalException,
117 "Cross-section biasing available only for gamma and electro nuclear reactions.");
118 }
119 if(aScale<100)
120 {
121 G4Exception("G4HadronicProcess", "001", JustWarning,
122 "Cross-section bias readjusted to be above safe limit. New value is 100");
123 aScaleFactor = 100.;
124 }
125 }
126
127 protected:
128
129 virtual void ResetNumberOfInteractionLengthLeft()
130 {
131 G4VProcess::theNumberOfInteractionLengthLeft =
132 -std::log( G4UniformRand() );
133 theInitialNumberOfInteractionLength =
134 G4VProcess::theNumberOfInteractionLengthLeft;
135 }
136
137 G4VParticleChange *GeneralPostStepDoIt( const G4Track &aTrack,
138 const G4Step &aStep );
139
140 void SetDispatch( G4HadronicProcess *value )
141 { dispatch=value; }
142
143 G4Element* ChooseAandZ(const G4DynamicParticle *aParticle,
144 const G4Material *aMaterial);
145
146 inline const G4EnergyRangeManager &GetEnergyRangeManager() const
147 { return theEnergyRangeManager; }
148
149 inline void SetEnergyRangeManager( const G4EnergyRangeManager &value )
150 { theEnergyRangeManager = value; }
151
152 inline G4HadronicInteraction *ChooseHadronicInteraction(
153 G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement )
154 {
155 G4EnergyRangeManager* ERMan = GetManagerPointer();
156 if(!ERMan->GetHadronicInteractionCounter())
157 G4cout<< "*G4HadronicProcess::ChooseHadronicInteraction: process = "
158 << GetProcessName() << ", nM="
159 << ERMan->GetHadronicInteractionCounter() << G4endl;
160 return ERMan->GetHadronicInteraction(kineticEnergy,aMaterial,anElement);
161 }
162
163 inline G4HadronicInteraction *GetHadronicInteraction()
164 { return theInteraction; }
165
166 public:
167 inline G4EnergyRangeManager *GetManagerPointer()
168 { return &theEnergyRangeManager; }
169 protected:
170
171 G4CrossSectionDataStore* GetCrossSectionDataStore()
172 {
173 return theCrossSectionDataStore;
174 }
175
176 G4double GetLastCrossSection() {return theLastCrossSection;}
177 private:
178
179 G4HadFinalState * DoIsotopeCounting(G4HadFinalState * aResult,
180 const G4Track & aTrack,
181 const G4Nucleus & aNucleus);
182
183 G4IsoResult * ExtractResidualNucleus(const G4Track & aTrack,
184 const G4Nucleus & aNucleus,
185 G4HadFinalState * aResult);
186
187 G4double GetTotalNumberOfInteractionLengthTraversed()
188 {
189 return theInitialNumberOfInteractionLength
190 -G4VProcess::theNumberOfInteractionLengthLeft;
191 }
192
193
194 void FillTotalResult(G4HadFinalState * aR, const G4Track & aT);
195
196 void SetCrossSectionDataStore(G4CrossSectionDataStore* aDataStore)
197 {
198 theCrossSectionDataStore = aDataStore;
199 }
200
201 G4double XBiasSurvivalProbability();
202 G4double XBiasSecondaryWeight();
203
204 private:
205
206 G4EnergyRangeManager theEnergyRangeManager;
207
208 G4HadronicInteraction *theInteraction;
209
210 G4CrossSectionDataStore* theCrossSectionDataStore;
211
212 G4Nucleus targetNucleus;
213
214 G4HadronicProcess *dispatch;
215
216// swiches for isotope production
217
218 static G4bool isoIsEnabled; // true or false; local swich overrides
219 G4int isoIsOnAnyway; // true(1), false(-1) or default(0)
220
221 G4IsoParticleChange theIsoPC;
222 std::vector<G4VIsotopeProduction *> theProductionModels;
223
224 std::vector<G4VLeadingParticleBiasing *> theBias;
225
226 static G4IsoParticleChange* theIsoResult;
227 static G4IsoParticleChange* theOldIsoResult;
228
229 G4ParticleChange* theTotalResult;
230
231 G4double theInitialNumberOfInteractionLength;
232
233 G4double aScaleFactor;
234 G4bool xBiasOn;
235 G4double theLastCrossSection;
236
237 G4int ModelingState;
238 };
239
240#endif
241
Note: See TracBrowser for help on using the repository browser.