source: trunk/source/processes/hadronic/models/management/include/G4HadronicInteraction.hh@ 1330

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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