Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4HadronicInteraction.cc
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// Hadronic Interaction base class
28// original by H.P. Wellisch
29// modified by J.L. Chuma, TRIUMF, 21-Mar-1997
30// Last modified: 04-Apr-1997
31// reimplemented 1.11.2003 JPW.
32// 23-Jan-2009 V.Ivanchenko move constructor and destructor to the body
33
34#include <iostream>
35
37#include "G4SystemOfUnits.hh"
40
42 verboseLevel(0), theMinEnergy(0.0),
43 isBlocked(false), recoilEnergyThreshold(0.0), theModelName(modelName),
44 epCheckLevels(DBL_MAX, DBL_MAX)
45{
48 registry->RegisterMe(this);
49}
50
52{
53 registry->RemoveMe(this);
54}
55
57{}
58
60{}
61
64{
65 return nullptr;
66}
67
71{
72 return 0.0;
73}
74
76 G4Nucleus&)
77{
78 return true;
79}
80
82 const G4Material *aMaterial, const G4Element *anElement ) const
83{
84 if(!IsBlocked()) { return theMinEnergy; }
85 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return DBL_MAX; }
86 if(!theMinEnergyListElements.empty()) {
87 for(auto const& elmlist : theMinEnergyListElements) {
88 if( anElement == elmlist.second )
89 { return elmlist.first; }
90 }
91 }
92 if(!theMinEnergyList.empty()) {
93 for(auto const & matlist : theMinEnergyList) {
94 if( aMaterial == matlist.second )
95 { return matlist.first; }
96 }
97 }
98 return theMinEnergy;
99}
100
102 const G4Element *anElement )
103{
104 Block();
105 if(!theMinEnergyListElements.empty()) {
106 for(auto & elmlist : theMinEnergyListElements) {
107 if( anElement == elmlist.second ) {
108 elmlist.first = anEnergy;
109 return;
110 }
111 }
112 }
113 theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
114}
115
117 const G4Material *aMaterial )
118{
119 Block();
120 if(!theMinEnergyList.empty()) {
121 for(auto & matlist : theMinEnergyList) {
122 if( aMaterial == matlist.second ) {
123 matlist.first = anEnergy;
124 return;
125 }
126 }
127 }
128 theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
129}
130
132 const G4Element *anElement ) const
133{
134 if(!IsBlocked()) { return theMaxEnergy; }
135 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return 0.0; }
136 if(!theMaxEnergyListElements.empty()) {
137 for(auto const& elmlist : theMaxEnergyListElements) {
138 if( anElement == elmlist.second )
139 { return elmlist.first; }
140 }
141 }
142 if(!theMaxEnergyList.empty()) {
143 for(auto const& matlist : theMaxEnergyList) {
144 if( aMaterial == matlist.second )
145 { return matlist.first; }
146 }
147 }
148 return theMaxEnergy;
149}
150
152 const G4Element *anElement )
153{
154 Block();
155 if(!theMaxEnergyListElements.empty()) {
156 for(auto & elmlist : theMaxEnergyListElements) {
157 if( anElement == elmlist.second ) {
158 elmlist.first = anEnergy;
159 return;
160 }
161 }
162 }
163 theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
164}
165
167{
168 Block();
169 if(!theMaxEnergyList.empty()) {
170 for(auto & matlist: theMaxEnergyList) {
171 if( aMaterial == matlist.second ) {
172 matlist.first = anEnergy;
173 return;
174 }
175 }
176 }
177 theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
178}
179
181{
182 Block();
183 theBlockedList.push_back(aMaterial);
184}
185
187{
188 Block();
189 theBlockedListElements.push_back(anElement);
190}
191
192
194{
195 for (auto const& mat : theBlockedList) {
196 if (aMaterial == mat) return true;
197 }
198 return false;
199}
200
201
203{
204 for (auto const& elm : theBlockedListElements) {
205 if (anElement == elm) return true;
206 }
207 return false;
208}
209
210const std::pair<G4double, G4double> G4HadronicInteraction::GetFatalEnergyCheckLevels() const
211{
212 // default level of Check
213 return std::pair<G4double, G4double>(2.*perCent, 1. * GeV);
214}
215
216std::pair<G4double, G4double>
218{
219 return epCheckLevels;
220}
221
222void G4HadronicInteraction::ModelDescription(std::ostream& outFile) const
223{
224 outFile << "The description for this model has not been written yet.\n";
225}
226
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void RegisterMe(G4HadronicInteraction *aModel)
void RemoveMe(G4HadronicInteraction *aModel)
static G4HadronicInteractionRegistry * Instance()
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void DeActivateFor(const G4Material *aMaterial)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
void SetMinEnergy(G4double anEnergy)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4HadronicInteraction(const G4String &modelName="HadronicModel")
virtual void ModelDescription(std::ostream &outFile) const
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SetMaxEnergy(const G4double anEnergy)
virtual G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
static G4HadronicParameters * Instance()
G4double GetMaxEnergy() const
#define DBL_MAX
Definition: templates.hh:62