Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4HadronicInteraction.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 //
27 // $Id: G4HadronicInteraction.hh 96490 2016-04-19 06:57:04Z gcosmo $
28 //
29 // Hadronic Interaction abstract base class
30 // This class is the base class for the model classes.
31 // It sorts out the energy-range for the models and provides
32 // class utilities.
33 // original by H.P. Wellisch
34 // Modified by J.L.Chuma, TRIUMF, 21-Mar-1997
35 // Last modified: 3-Apr-1997
36 // Added units to energy initialization: J.L. Chuma 04-Apr-97
37 // Modified by J.L.Chuma, 05-May-97 to Initialize theBlockedCounter
38 // Modified by J.L.Chuma, 08-Jul-97 to implement the Nucleus changes
39 // Adding a registry for memory management of hadronic models, HPW 22-Mar-99
40 // 23-Jan-2009 V.Ivanchenko move constructor and destructor to the body
41 // and reorder methods in the header
42 // 29-Jun-2009 V.Ivanchenko add SampleInvariantT method
43 // 29-Aug-2009 V.Ivanchenko moved G4ReactionDynamics to G4InelasticInteraction,
44 // add const pointers, and added recoilEnergyThreshold
45 // member and accesors
46 
47 // Class Description
48 // This is the base class for all hadronic interaction models in geant4.
49 // If you want to implement a new way of producing a final state, please,
50 // inherit from here.
51 // Class Description - End
52 
53 #ifndef G4HadronicInteraction_h
54 #define G4HadronicInteraction_h 1
55 
56 #include "G4HadFinalState.hh"
57 #include "G4Material.hh"
58 #include "G4Nucleus.hh"
59 #include "G4Track.hh"
60 #include "G4HadProjectile.hh"
61 
63 
65 {
66 public: // With description
67 
68  explicit G4HadronicInteraction(const G4String& modelName = "HadronicModel");
69 
70  virtual ~G4HadronicInteraction();
71 
72  virtual G4HadFinalState *ApplyYourself(const G4HadProjectile &aTrack,
73  G4Nucleus & targetNucleus ) = 0;
74  // The interface to implement for final state production code.
75 
77  G4double plab,
78  G4int Z, G4int A);
79  // The interface to implement sampling of scattering or change exchange
80 
81  virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
82  G4Nucleus & targetNucleus);
83 
84  inline G4double GetMinEnergy() const
85  { return theMinEnergy; }
86 
87  G4double GetMinEnergy( const G4Material *aMaterial,
88  const G4Element *anElement ) const;
89 
90  inline void SetMinEnergy( G4double anEnergy )
91  { theMinEnergy = anEnergy; }
92 
93  void SetMinEnergy( G4double anEnergy, const G4Element *anElement );
94 
95  void SetMinEnergy( G4double anEnergy, const G4Material *aMaterial );
96 
97  inline G4double GetMaxEnergy() const
98  { return theMaxEnergy; }
99 
100  G4double GetMaxEnergy( const G4Material *aMaterial,
101  const G4Element *anElement ) const;
102 
103  inline void SetMaxEnergy( const G4double anEnergy )
104  { theMaxEnergy = anEnergy; }
105 
106  void SetMaxEnergy( G4double anEnergy, const G4Element *anElement );
107 
108  void SetMaxEnergy( G4double anEnergy, const G4Material *aMaterial );
109 
110  inline G4int GetVerboseLevel() const
111  { return verboseLevel; }
112 
113  inline void SetVerboseLevel( G4int value )
114  { verboseLevel = value; }
115 
116  inline const G4String& GetModelName() const
117  { return theModelName; }
118 
119  void DeActivateFor(const G4Material* aMaterial);
120 
121  inline void ActivateFor( const G4Material *aMaterial )
122  {
123  Block();
124  SetMaxEnergy(GetMaxEnergy(), aMaterial);
125  SetMinEnergy(GetMinEnergy(), aMaterial);
126  }
127 
128  void DeActivateFor( const G4Element *anElement );
129  inline void ActivateFor( const G4Element *anElement )
130  {
131  Block();
132  SetMaxEnergy(GetMaxEnergy(), anElement);
133  SetMinEnergy(GetMinEnergy(), anElement);
134  }
135 
136  G4bool IsBlocked( const G4Material *aMaterial ) const;
137  G4bool IsBlocked( const G4Element *anElement) const;
138 
140  { recoilEnergyThreshold = val; }
141 
143  { return recoilEnergyThreshold;}
144 
145  virtual const std::pair<G4double, G4double> GetFatalEnergyCheckLevels() const;
146 
147  virtual std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const;
148 
149  inline void
150  SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
151  { epCheckLevels.first = relativeLevel;
152  epCheckLevels.second = absoluteLevel; }
153 
154  virtual void ModelDescription(std::ostream& outFile) const ; //=0;
155 
156  // Initialisation before run
157  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
158  virtual void InitialiseModel();
159 
160 protected:
161 
162  inline void SetModelName(const G4String& nam)
163  { theModelName = nam; }
164 
165  inline G4bool IsBlocked() const { return isBlocked;}
166  inline void Block() { isBlocked = true; }
167 
169  // the G4HadFinalState object which is modified and returned
170  // by address by the ApplyYourself method,
171  // (instead of aParticleChange as found in G4VProcess)
172 
174  // control flag for output messages
175  // 0: silent
176  // 1: warning messages
177  // 2: more
178  // (instead of verboseLevel as found in G4VProcess)
179 
180  // these two have global validity energy range
183 
185 
186 private:
187 
189  const G4HadronicInteraction& operator=(const G4HadronicInteraction &right) = delete;
190  G4bool operator==(const G4HadronicInteraction &right ) const = delete;
191  G4bool operator!=(const G4HadronicInteraction &right ) const = delete;
192 
194 
195  G4double recoilEnergyThreshold;
196 
197  G4String theModelName;
198 
199  std::pair<G4double, G4double> epCheckLevels;
200 
201  std::vector<std::pair<G4double, const G4Material *> > theMinEnergyList;
202  std::vector<std::pair<G4double, const G4Material *> > theMaxEnergyList;
203  std::vector<std::pair<G4double, const G4Element *> > theMinEnergyListElements;
204  std::vector<std::pair<G4double, const G4Element *> > theMaxEnergyListElements;
205  std::vector<const G4Material *> theBlockedList;
206  std::vector<const G4Element *> theBlockedListElements;
207 };
208 
209 #endif
G4double GetMinEnergy() const
virtual G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
G4double GetMaxEnergy() const
const char * p
Definition: xmltok.h:285
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
virtual void ModelDescription(std::ostream &outFile) const
const G4String & GetModelName() const
int G4int
Definition: G4Types.hh:78
G4double GetRecoilEnergyThreshold() const
void SetMinEnergy(G4double anEnergy)
double A(double temperature)
const XML_Char int const XML_Char * value
Definition: expat.h:331
bool G4bool
Definition: G4Types.hh:79
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
void ActivateFor(const G4Material *aMaterial)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetModelName(const G4String &nam)
void SetMaxEnergy(const G4double anEnergy)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void SetRecoilEnergyThreshold(G4double val)
double G4double
Definition: G4Types.hh:76
void SetVerboseLevel(G4int value)
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
void ActivateFor(const G4Element *anElement)
void DeActivateFor(const G4Material *aMaterial)