Geant4  10.01.p02
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 83773 2014-09-15 07:19:09Z 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 {
64 public: // With description
65 
66  G4HadronicInteraction(const G4String& modelName = "HadronicModel");
67 
68  virtual ~G4HadronicInteraction();
69 
70  virtual G4HadFinalState *ApplyYourself(const G4HadProjectile &aTrack,
71  G4Nucleus & targetNucleus ) = 0;
72  // The interface to implement for final state production code.
73 
75  G4double plab,
76  G4int Z, G4int A);
77  // The interface to implement sampling of scattering or change exchange
78 
79  virtual G4bool IsApplicable(const G4HadProjectile &/*aTrack*/,
80  G4Nucleus & /*targetNucleus*/)
81  { return true;}
82 
83  inline G4double GetMinEnergy() const
84  { return theMinEnergy; }
85 
86  G4double GetMinEnergy( const G4Material *aMaterial,
87  const G4Element *anElement ) const;
88 
89  inline void SetMinEnergy( G4double anEnergy )
90  { theMinEnergy = anEnergy; }
91 
92  void SetMinEnergy( G4double anEnergy, const G4Element *anElement );
93 
94  void SetMinEnergy( G4double anEnergy, const G4Material *aMaterial );
95 
96  inline G4double GetMaxEnergy() const
97  { return theMaxEnergy; }
98 
99  G4double GetMaxEnergy( const G4Material *aMaterial,
100  const G4Element *anElement ) const;
101 
102  inline void SetMaxEnergy( const G4double anEnergy )
103  { theMaxEnergy = anEnergy; }
104 
105  void SetMaxEnergy( G4double anEnergy, const G4Element *anElement );
106 
107  void SetMaxEnergy( G4double anEnergy, const G4Material *aMaterial );
108 
109  inline const G4HadronicInteraction* GetMyPointer() const
110  { return this; }
111 
112  virtual G4int GetVerboseLevel() const
113  { return verboseLevel; }
114 
115  virtual void SetVerboseLevel( G4int value )
116  { verboseLevel = value; }
117 
118  inline const G4String& GetModelName() const
119  { return theModelName; }
120 
121  void DeActivateFor(const G4Material* aMaterial);
122 
123  inline void ActivateFor( const G4Material *aMaterial )
124  {
125  Block();
126  SetMaxEnergy(GetMaxEnergy(), aMaterial);
127  SetMinEnergy(GetMinEnergy(), aMaterial);
128  }
129 
130  void DeActivateFor( const G4Element *anElement );
131  inline void ActivateFor( const G4Element *anElement )
132  {
133  Block();
134  SetMaxEnergy(GetMaxEnergy(), anElement);
135  SetMinEnergy(GetMinEnergy(), anElement);
136  }
137 
138  G4bool IsBlocked( const G4Material *aMaterial ) const;
139  G4bool IsBlocked( const G4Element *anElement) const;
140 
142  { recoilEnergyThreshold = val; }
143 
145  { return recoilEnergyThreshold;}
146 
148  { return ( this == (G4HadronicInteraction *) &right ); }
149 
151  { return ( this != (G4HadronicInteraction *) &right ); }
152 
153  virtual const std::pair<G4double, G4double> GetFatalEnergyCheckLevels() const;
154 
155  virtual std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const;
156 
157  inline void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
158  { epCheckLevels.first = relativeLevel;
159  epCheckLevels.second = absoluteLevel; }
160 
161  virtual void ModelDescription(std::ostream& outFile) const ; //=0;
162 
163  // Initialisation before run
164  virtual void BuildPhysicsTable(const G4ParticleDefinition&){;};
165 
166 private:
167 
170 
171 protected:
172 
173  inline void SetModelName(const G4String& nam)
174  { theModelName = nam; }
175 
176  inline G4bool IsBlocked() const { return isBlocked;}
177  inline void Block() { isBlocked = true; }
178 
180  // the G4HadFinalState object which is modified and returned
181  // by address by the ApplyYourself method,
182  // (instead of aParticleChange as found in G4VProcess)
183 
185  // control flag for output messages
186  // 0: silent
187  // 1: warning messages
188  // 2: more
189  // (instead of verboseLevel as found in G4VProcess)
190 
191  // these two have global validity energy range
194 
196 
197 
198 private:
199 
201 
203 
204  std::pair<G4double, G4double> epCheckLevels;
205 
206  std::vector<std::pair<G4double, const G4Material *> > theMinEnergyList;
207  std::vector<std::pair<G4double, const G4Material *> > theMaxEnergyList;
208  std::vector<std::pair<G4double, const G4Element *> > theMinEnergyListElements;
209  std::vector<std::pair<G4double, const G4Element *> > theMaxEnergyListElements;
210  std::vector<const G4Material *> theBlockedList;
211  std::vector<const G4Element *> theBlockedListElements;
212 };
213 
214 #endif
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements
G4double GetMinEnergy() const
std::vector< const G4Element * > theBlockedListElements
const G4HadronicInteraction * GetMyPointer() const
virtual G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
virtual G4bool IsApplicable(const G4HadProjectile &, G4Nucleus &)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetMaxEnergy() const
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
virtual void ModelDescription(std::ostream &outFile) const
G4bool operator!=(const G4HadronicInteraction &right) const
const G4String & GetModelName() const
int G4int
Definition: G4Types.hh:78
G4double GetRecoilEnergyThreshold() const
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
void SetMinEnergy(G4double anEnergy)
virtual G4int GetVerboseLevel() const
bool G4bool
Definition: G4Types.hh:79
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
void ActivateFor(const G4Material *aMaterial)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static const G4double A[nN]
std::pair< G4double, G4double > epCheckLevels
void SetModelName(const G4String &nam)
const G4HadronicInteraction & operator=(const G4HadronicInteraction &right)
virtual void SetVerboseLevel(G4int value)
void SetMaxEnergy(const G4double anEnergy)
std::vector< const G4Material * > theBlockedList
G4bool operator==(const G4HadronicInteraction &right) const
void SetRecoilEnergyThreshold(G4double val)
double G4double
Definition: G4Types.hh:76
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)