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

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

update processes

File size: 7.3 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
59class G4HadronicProcess : public G4VDiscreteProcess
60{
61public:
62   
63  G4HadronicProcess( const G4String &processName = "Hadronic", 
64                     G4ProcessType   aType = fHadronic );   
65
66  virtual ~G4HadronicProcess();
67
68  // register generator of secondaries
69  void RegisterMe( G4HadronicInteraction *a );
70
71  // get cross section per element
72  virtual 
73  G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle, 
74                                      const G4Element *anElement, 
75                                      G4double aTemp );
76
77  // generic PostStepDoIt recommended for all derived classes
78  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
79                                          const G4Step& aStep);
80
81  // initialisation of physics tables and G4HadronicProcessStore
82  virtual void PreparePhysicsTable(const G4ParticleDefinition&);
83
84  // build physics tables and print out the configuration of the process
85  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
86
87  // dump physics tables
88  inline void DumpPhysicsTable(const G4ParticleDefinition& p)
89  { theCrossSectionDataStore->DumpPhysicsTable(p); }
90
91  // add cross section data set
92  inline void AddDataSet(G4VCrossSectionDataSet * aDataSet)
93  { theCrossSectionDataStore->AddDataSet(aDataSet);}
94
95  // access to the manager
96  inline G4EnergyRangeManager *GetManagerPointer()
97  { return &theEnergyRangeManager; }
98         
99  // get inverse cross section per volume
100  G4double GetMeanFreePath(const G4Track &aTrack, G4double, 
101                           G4ForceCondition *);
102
103protected:   
104
105  // reset number of interaction length and save 
106  virtual void ResetNumberOfInteractionLengthLeft()
107  { G4VProcess::ResetNumberOfInteractionLengthLeft(); 
108    theInitialNumberOfInteractionLength = 
109      G4VProcess::theNumberOfInteractionLengthLeft;
110  }
111
112  // generic method to choose secondary generator
113  // recommended for all derived classes
114  inline G4HadronicInteraction *ChooseHadronicInteraction(
115      G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement )
116  { return theEnergyRangeManager.GetHadronicInteraction(kineticEnergy,
117                                                        aMaterial,anElement);
118  }
119
120public:
121
122  // Methods for isotope production   
123  static void EnableIsotopeProductionGlobally();
124  static void DisableIsotopeProductionGlobally();
125   
126  void EnableIsotopeCounting()  {isoIsOnAnyway = 1;}
127  void DisableIsotopeCounting() {isoIsOnAnyway = -1;}
128   
129  void RegisterIsotopeProductionModel(G4VIsotopeProduction * aModel)
130  { theProductionModels.push_back(aModel); }
131
132  static G4IsoParticleChange * GetIsotopeProductionInfo();
133
134  void BiasCrossSectionByFactor(G4double aScale);
135   
136protected:
137           
138  // obsolete method will be removed
139  inline const G4EnergyRangeManager &GetEnergyRangeManager() const
140  { return theEnergyRangeManager; }
141   
142  // obsolete method will be removed
143  inline void SetEnergyRangeManager( const G4EnergyRangeManager &value )
144  { theEnergyRangeManager = value; }
145
146  // access to the chosen generator
147  inline G4HadronicInteraction *GetHadronicInteraction()
148  { return theInteraction; }
149   
150  // access to the cross section data store
151  inline G4CrossSectionDataStore* GetCrossSectionDataStore()
152  { return theCrossSectionDataStore; }
153   
154  // access to the cross section data set
155  inline G4double GetLastCrossSection() 
156  { return theLastCrossSection; }
157
158private:
159   
160  void FillTotalResult(G4HadFinalState * aR, const G4Track & aT);
161
162  G4HadFinalState * DoIsotopeCounting(G4HadFinalState * aResult,
163                                      const G4Track & aTrack,
164                                      const G4Nucleus & aNucleus);
165                                         
166  G4IsoResult * ExtractResidualNucleus(const G4Track & aTrack,
167                                       const G4Nucleus & aNucleus,
168                                       G4HadFinalState * aResult);
169
170  inline G4double GetTotalNumberOfInteractionLengthTraversed()
171  { return theInitialNumberOfInteractionLength
172      -G4VProcess::theNumberOfInteractionLengthLeft;
173  }
174       
175  //  inline void SetCrossSectionDataStore(G4CrossSectionDataStore* aDataStore)
176  //  { theCrossSectionDataStore = aDataStore; }
177   
178  G4double XBiasSurvivalProbability();
179  G4double XBiasSecondaryWeight();
180   
181private:
182   
183  G4EnergyRangeManager theEnergyRangeManager;
184   
185  G4HadronicInteraction *theInteraction;
186
187  G4CrossSectionDataStore* theCrossSectionDataStore;
188 
189  G4Nucleus targetNucleus;
190   
191  G4HadronicProcess *dispatch;
192
193  bool G4HadronicProcess_debug_flag;
194
195  // swiches for isotope production   
196  static G4bool isoIsEnabled; // true or false; local swich overrides
197  G4int isoIsOnAnyway; // true(1), false(-1) or default(0)
198   
199  G4IsoParticleChange theIsoPC;
200  std::vector<G4VIsotopeProduction *> theProductionModels;
201   
202  std::vector<G4VLeadingParticleBiasing *> theBias;
203
204  static G4IsoParticleChange* theIsoResult;
205  static G4IsoParticleChange* theOldIsoResult;
206   
207  G4ParticleChange* theTotalResult; 
208   
209  G4double theInitialNumberOfInteractionLength;   
210
211  G4double aScaleFactor;
212  G4bool xBiasOn;
213  G4double theLastCrossSection;
214
215  G4int ModelingState;
216};
217 
218#endif
219 
Note: See TracBrowser for help on using the repository browser.