Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id: G4HadronicInteraction.cc 96490 2016-04-19 06:57:04Z gcosmo $
27 //
28 // Hadronic Interaction base class
29 // original by H.P. Wellisch
30 // modified by J.L. Chuma, TRIUMF, 21-Mar-1997
31 // Last modified: 04-Apr-1997
32 // reimplemented 1.11.2003 JPW.
33 // 23-Jan-2009 V.Ivanchenko move constructor and destructor to the body
34 
35 #include <iostream>
36 
37 #include "G4HadronicInteraction.hh"
38 #include "G4SystemOfUnits.hh"
40 #include "G4HadronicException.hh"
41 
43  verboseLevel(0), theMinEnergy(0.0), theMaxEnergy(25.0*GeV),
44  isBlocked(false), recoilEnergyThreshold(0.0), theModelName(modelName),
45  epCheckLevels(DBL_MAX, DBL_MAX)
46 {
48  registry->RegisterMe(this);
49 }
50 
52 {
53  registry->RemoveMe(this);
54 }
55 
57 {}
58 
60 {}
61 
62 G4double
65 {
66  return 0.0;
67 }
68 
70  G4Nucleus&)
71 {
72  return true;
73 }
74 
75 
77  const G4Material *aMaterial, const G4Element *anElement ) const
78 {
79  if( IsBlocked(aMaterial) ) { return 0.0; }
80  if( IsBlocked(anElement) ) { return 0.0; }
81  size_t length = theMinEnergyListElements.size();
82  if(0 < length) {
83  for(size_t i=0; i<length; ++i ) {
84  if( anElement == theMinEnergyListElements[i].second )
85  { return theMinEnergyListElements[i].first; }
86  }
87  }
88  length = theMinEnergyList.size();
89  if(0 < length) {
90  for(size_t i=0; i<length; ++i ) {
91  if( aMaterial == theMinEnergyList[i].second )
92  { return theMinEnergyList[i].first; }
93  }
94  }
95  if(IsBlocked()) { return 0.0; }
96  if( verboseLevel > 1 ) {
97  G4cout << "*** Warning from HadronicInteraction::GetMinEnergy" << G4endl
98  << " material " << aMaterial->GetName()
99  << " not found in min energy List" << G4endl;
100  }
101  return theMinEnergy;
102 }
103 
105  const G4Element *anElement )
106 {
107  if( IsBlocked(anElement) ) {
108  G4cout << "*** Warning from HadronicInteraction::SetMinEnergy" << G4endl
109  << " The model is not active for the Element "
110  << anElement->GetName() << "." << G4endl;
111  }
112  size_t length = theMinEnergyListElements.size();
113  if(0 < length) {
114  for(size_t i=0; i<length; ++i ) {
115  if( anElement == theMinEnergyListElements[i].second )
116  {
117  theMinEnergyListElements[i].first = anEnergy;
118  return;
119  }
120  }
121  }
122  theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
123 }
124 
126  const G4Material *aMaterial )
127 {
128  if( IsBlocked(aMaterial) ) {
129  G4cout << "*** Warning from HadronicInteraction::SetMinEnergy" << G4endl
130  << " The model is not active for the Material "
131  << aMaterial->GetName() << "." << G4endl;
132  }
133  size_t length = theMinEnergyList.size();
134  if(0 < length) {
135  for(size_t i=0; i<length; ++i ) {
136  if( aMaterial == theMinEnergyList[i].second )
137  {
138  theMinEnergyList[i].first = anEnergy;
139  return;
140  }
141  }
142  }
143  theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
144 }
145 
147  const G4Element *anElement ) const
148 {
149  if( IsBlocked(aMaterial) ) { return 0.0; }
150  if( IsBlocked(anElement) ) { return 0.0; }
151  size_t length = theMaxEnergyListElements.size();
152  if(0 < length) {
153  for(size_t i=0; i<length; ++i ) {
154  if( anElement == theMaxEnergyListElements[i].second )
155  { return theMaxEnergyListElements[i].first; }
156  }
157  }
158  length = theMaxEnergyList.size();
159  if(0 < length) {
160  for(size_t i=0; i<length; ++i ) {
161  if( aMaterial == theMaxEnergyList[i].second )
162  { return theMaxEnergyList[i].first; }
163  }
164  }
165  if(IsBlocked()) { return 0.0; }
166  if( verboseLevel > 1 ) {
167  G4cout << "*** Warning from HadronicInteraction::GetMaxEnergy" << G4endl
168  << " material " << aMaterial->GetName()
169  << " not found in min energy List" << G4endl;
170  }
171  return theMaxEnergy;
172 }
173 
175  const G4Element *anElement )
176 {
177  if( IsBlocked(anElement) ) {
178  G4cout << "*** Warning from HadronicInteraction::SetMaxEnergy" << G4endl
179  << "Warning: The model is not active for the Element "
180  << anElement->GetName() << "." << G4endl;
181  }
182  size_t length = theMaxEnergyListElements.size();
183  if(0 < length) {
184  for(size_t i=0; i<length; ++i ) {
185  if( anElement == theMaxEnergyListElements[i].second )
186  {
187  theMaxEnergyListElements[i].first = anEnergy;
188  return;
189  }
190  }
191  }
192  theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
193 }
194 
196  const G4Material *aMaterial )
197 {
198  if( IsBlocked(aMaterial) ) {
199  G4cout << "*** Warning from HadronicInteraction::SetMaxEnergy" << G4endl
200  << "Warning: The model is not active for the Material "
201  << aMaterial->GetName() << "." << G4endl;
202  }
203  size_t length = theMaxEnergyList.size();
204  if(0 < length) {
205  for(size_t i=0; i<length; ++i ) {
206  if( aMaterial == theMaxEnergyList[i].second )
207  {
208  theMaxEnergyList[i].first = anEnergy;
209  return;
210  }
211  }
212  }
213  theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
214 }
215 
217 {
218  theBlockedList.push_back(aMaterial);
219 }
220 
222 {
223  theBlockedListElements.push_back(anElement);
224 }
225 
226 
228 {
229  for (size_t i=0; i<theBlockedList.size(); ++i) {
230  if (aMaterial == theBlockedList[i]) return true;
231  }
232  return false;
233 }
234 
235 
237 {
238  for (size_t i=0; i<theBlockedListElements.size(); ++i) {
239  if (anElement == theBlockedListElements[i]) return true;
240  }
241  return false;
242 }
243 
244 const std::pair<G4double, G4double> G4HadronicInteraction::GetFatalEnergyCheckLevels() const
245 {
246  // default level of Check
247  return std::pair<G4double, G4double>(2.*perCent, 1. * GeV);
248 }
249 
250 std::pair<G4double, G4double>
252 {
253  return epCheckLevels;
254 }
255 
256 
257 void G4HadronicInteraction::ModelDescription(std::ostream& outFile) const
258 {
259  outFile << "The description for this model has not been written yet.\n";
260 }
261 
262 /*
263 G4HadronicInteraction::G4HadronicInteraction(const G4HadronicInteraction &right )
264 {
265  *this = right;
266 }
267 
268 const G4HadronicInteraction&
269 G4HadronicInteraction::operator=(const G4HadronicInteraction &right )
270 {
271  G4String text = "unintended use of G4HadronicInteraction::operator=";
272  throw G4HadronicException(__FILE__, __LINE__, text);
273  return right;
274 }
275  */
276 /* end of file */
277 
G4double GetMinEnergy() const
G4bool isBlocked(ParticleList const &p, Nucleus const *const n)
Check Pauli blocking.
virtual G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
G4double GetMaxEnergy() const
static constexpr double perCent
Definition: G4SIunits.hh:332
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
const G4String & GetName() const
Definition: G4Material.hh:178
virtual void ModelDescription(std::ostream &outFile) const
void RemoveMe(G4HadronicInteraction *aModel)
static constexpr double second
Definition: G4SIunits.hh:157
void RegisterMe(G4HadronicInteraction *aModel)
int G4int
Definition: G4Types.hh:78
void SetMinEnergy(G4double anEnergy)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4HadronicInteractionRegistry * Instance()
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
#define DBL_MAX
Definition: templates.hh:83
void DeActivateFor(const G4Material *aMaterial)