Geant4  10.02.p01
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 88839 2015-03-12 10:30:07Z 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 
56 G4double
59 {
60  return 0.0;
61 }
62 
64  const G4Material *aMaterial, const G4Element *anElement ) const
65 {
66  if( IsBlocked(aMaterial) ) { return 0.0; }
67  if( IsBlocked(anElement) ) { return 0.0; }
68  size_t length = theMinEnergyListElements.size();
69  if(0 < length) {
70  for(size_t i=0; i<length; ++i ) {
71  if( anElement == theMinEnergyListElements[i].second )
72  { return theMinEnergyListElements[i].first; }
73  }
74  }
75  length = theMinEnergyList.size();
76  if(0 < length) {
77  for(size_t i=0; i<length; ++i ) {
78  if( aMaterial == theMinEnergyList[i].second )
79  { return theMinEnergyList[i].first; }
80  }
81  }
82  if(IsBlocked()) { return 0.0; }
83  if( verboseLevel > 1 ) {
84  G4cout << "*** Warning from HadronicInteraction::GetMinEnergy" << G4endl
85  << " material " << aMaterial->GetName()
86  << " not found in min energy List" << G4endl;
87  }
88  return theMinEnergy;
89 }
90 
92  const G4Element *anElement )
93 {
94  if( IsBlocked(anElement) ) {
95  G4cout << "*** Warning from HadronicInteraction::SetMinEnergy" << G4endl
96  << " The model is not active for the Element "
97  << anElement->GetName() << "." << G4endl;
98  }
99  size_t length = theMinEnergyListElements.size();
100  if(0 < length) {
101  for(size_t i=0; i<length; ++i ) {
102  if( anElement == theMinEnergyListElements[i].second )
103  {
104  theMinEnergyListElements[i].first = anEnergy;
105  return;
106  }
107  }
108  }
109  theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
110 }
111 
113  const G4Material *aMaterial )
114 {
115  if( IsBlocked(aMaterial) ) {
116  G4cout << "*** Warning from HadronicInteraction::SetMinEnergy" << G4endl
117  << " The model is not active for the Material "
118  << aMaterial->GetName() << "." << G4endl;
119  }
120  size_t length = theMinEnergyList.size();
121  if(0 < length) {
122  for(size_t i=0; i<length; ++i ) {
123  if( aMaterial == theMinEnergyList[i].second )
124  {
125  theMinEnergyList[i].first = anEnergy;
126  return;
127  }
128  }
129  }
130  theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
131 }
132 
134  const G4Element *anElement ) const
135 {
136  if( IsBlocked(aMaterial) ) { return 0.0; }
137  if( IsBlocked(anElement) ) { return 0.0; }
138  size_t length = theMaxEnergyListElements.size();
139  if(0 < length) {
140  for(size_t i=0; i<length; ++i ) {
141  if( anElement == theMaxEnergyListElements[i].second )
142  { return theMaxEnergyListElements[i].first; }
143  }
144  }
145  length = theMaxEnergyList.size();
146  if(0 < length) {
147  for(size_t i=0; i<length; ++i ) {
148  if( aMaterial == theMaxEnergyList[i].second )
149  { return theMaxEnergyList[i].first; }
150  }
151  }
152  if(IsBlocked()) { return 0.0; }
153  if( verboseLevel > 1 ) {
154  G4cout << "*** Warning from HadronicInteraction::GetMaxEnergy" << G4endl
155  << " material " << aMaterial->GetName()
156  << " not found in min energy List" << G4endl;
157  }
158  return theMaxEnergy;
159 }
160 
162  const G4Element *anElement )
163 {
164  if( IsBlocked(anElement) ) {
165  G4cout << "*** Warning from HadronicInteraction::SetMaxEnergy" << G4endl
166  << "Warning: The model is not active for the Element "
167  << anElement->GetName() << "." << G4endl;
168  }
169  size_t length = theMaxEnergyListElements.size();
170  if(0 < length) {
171  for(size_t i=0; i<length; ++i ) {
172  if( anElement == theMaxEnergyListElements[i].second )
173  {
174  theMaxEnergyListElements[i].first = anEnergy;
175  return;
176  }
177  }
178  }
179  theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
180 }
181 
183  const G4Material *aMaterial )
184 {
185  if( IsBlocked(aMaterial) ) {
186  G4cout << "*** Warning from HadronicInteraction::SetMaxEnergy" << G4endl
187  << "Warning: The model is not active for the Material "
188  << aMaterial->GetName() << "." << G4endl;
189  }
190  size_t length = theMaxEnergyList.size();
191  if(0 < length) {
192  for(size_t i=0; i<length; ++i ) {
193  if( aMaterial == theMaxEnergyList[i].second )
194  {
195  theMaxEnergyList[i].first = anEnergy;
196  return;
197  }
198  }
199  }
200  theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
201 }
202 
204 {
205  theBlockedList.push_back(aMaterial);
206 }
207 
209 {
210  theBlockedListElements.push_back(anElement);
211 }
212 
213 
215 {
216  for (size_t i=0; i<theBlockedList.size(); ++i) {
217  if (aMaterial == theBlockedList[i]) return true;
218  }
219  return false;
220 }
221 
222 
224 {
225  for (size_t i=0; i<theBlockedListElements.size(); ++i) {
226  if (anElement == theBlockedListElements[i]) return true;
227  }
228  return false;
229 }
230 
231 const std::pair<G4double, G4double> G4HadronicInteraction::GetFatalEnergyCheckLevels() const
232 {
233  // default level of Check
234  return std::pair<G4double, G4double>(2.*perCent, 1. * GeV);
235 }
236 
237 std::pair<G4double, G4double>
239 {
240  return epCheckLevels;
241 }
242 
243 
244 void G4HadronicInteraction::ModelDescription(std::ostream& outFile) const
245 {
246  outFile << "The description for this model has not been written yet.\n";
247 }
248 
249 /*
250 G4HadronicInteraction::G4HadronicInteraction(const G4HadronicInteraction &right )
251 {
252  *this = right;
253 }
254 
255 const G4HadronicInteraction&
256 G4HadronicInteraction::operator=(const G4HadronicInteraction &right )
257 {
258  G4String text = "unintended use of G4HadronicInteraction::operator=";
259  throw G4HadronicException(__FILE__, __LINE__, text);
260  return right;
261 }
262  */
263 /* end of file */
264 
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements
G4double GetMinEnergy() const
G4bool isBlocked(ParticleList const &p, Nucleus const *const n)
Check Pauli blocking.
std::vector< const G4Element * > theBlockedListElements
virtual G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
G4double GetMaxEnergy() const
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)
G4HadronicInteractionRegistry * registry
void RegisterMe(G4HadronicInteraction *aModel)
int G4int
Definition: G4Types.hh:78
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
void SetMinEnergy(G4double anEnergy)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
static const double second
Definition: G4SIunits.hh:156
static const double GeV
Definition: G4SIunits.hh:214
static const double perCent
Definition: G4SIunits.hh:329
G4HadronicInteraction(const G4String &modelName="HadronicModel")
std::pair< G4double, G4double > epCheckLevels
static G4HadronicInteractionRegistry * Instance()
void SetMaxEnergy(const G4double anEnergy)
std::vector< const G4Material * > theBlockedList
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
#define DBL_MAX
Definition: templates.hh:83
void DeActivateFor(const G4Material *aMaterial)