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