Geant4  10.02.p02
G4EnergyRangeManager.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: G4EnergyRangeManager.cc 83772 2014-09-15 07:18:08Z gcosmo $
27 //
28  // Hadronic Process: Energy Range Manager
29  // original by H.P. Wellisch
30  // modified by J.L. Chuma, TRIUMF, 22-Nov-1996
31  // Last modified: 24-Mar-1997
32  // fix in the counter-hndling: H.P. Wellisch 04-Apr-97
33  // throw an exception if no model found: J.L. Chuma 04-Apr-97
34 
35 #include "G4EnergyRangeManager.hh"
36 #include "Randomize.hh"
37 #include "G4HadronicException.hh"
38 
40  : theHadronicInteractionCounter(0)
41 {}
42 
44 {}
45 
47 {
48  if (this != &right) {
51  }
52 }
53 
56 {
57  if (this != &right) {
60  }
61  return *this;
62 }
63 
65 {
66  if(!a) { return; }
68  for(G4int i=0; i<theHadronicInteractionCounter; ++i) {
69  if(a == theHadronicInteraction[i]) { return; }
70  }
71  }
72  theHadronicInteraction.push_back(a);
74 }
75 
76 
79  G4Nucleus & aTargetNucleus,
80  const G4Material* aMaterial,
81  const G4Element* anElement) const
82 {
84  throw G4HadronicException(__FILE__, __LINE__,
85  "GetHadronicInteraction: NO MODELS STORED");
86  }
87 
88  G4double kineticEnergy = aHadProjectile.GetKineticEnergy();
89  // For ions, get kinetic energy per nucleon
90  if ( aHadProjectile.GetDefinition()->GetBaryonNumber() > 1.5 ) {
91  kineticEnergy /= aHadProjectile.GetDefinition()->GetBaryonNumber();
92  }
93 
94  G4int cou = 0, memory = 0, memor2 = 0;
95  G4double emi1 = 0.0, ema1 = 0.0, emi2 = 0.0, ema2 = 0.0;
96 
97  for (G4int i = 0; i<theHadronicInteractionCounter; ++i) {
98  if ( theHadronicInteraction[i]->IsApplicable( aHadProjectile, aTargetNucleus ) ) {
99  G4double low = theHadronicInteraction[i]->GetMinEnergy( aMaterial, anElement );
100  // Work-around for particles with 0 kinetic energy, which still
101  // require a model to return a ParticleChange
102  //if (low == 0.) low = -DBL_MIN;
103  G4double high = theHadronicInteraction[i]->GetMaxEnergy( aMaterial, anElement );
104  if (low <= kineticEnergy && high > kineticEnergy) {
105  ++cou;
106  emi2 = emi1;
107  ema2 = ema1;
108  emi1 = low;
109  ema1 = high;
110  memor2 = memory;
111  memory = i;
112  }
113  }
114  }
115 
116  G4int mem = -1;
117  G4double rand;
118  switch (cou) {
119  case 0:
120  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="
121  <<theHadronicInteractionCounter<<", Ek="
122  <<kineticEnergy<<", Material = "<<aMaterial->GetName()
123  <<", Element = "
124  <<anElement->GetName()<<G4endl;
125  for( G4int j=0; j<theHadronicInteractionCounter; ++j)
126  {
128  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
129  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
130  }
131  throw G4HadronicException(__FILE__, __LINE__,
132  "GetHadronicInteraction: No Model found");
133  return 0;
134  case 1:
135  mem = memory;
136  break;
137  case 2:
138  if( (emi2<=emi1 && ema2>=ema1) || (emi2>=emi1 && ema2<=ema1) )
139  {
140  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="
141  <<theHadronicInteractionCounter<<", Ek="
142  <<kineticEnergy<<", Material = "<<aMaterial->GetName()
143  <<", Element = "
144  <<anElement->GetName()<<G4endl;
145  for( G4int j=0; j<theHadronicInteractionCounter; ++j)
146  {
148  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
149  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
150  }
151  throw G4HadronicException(__FILE__, __LINE__,
152  "GetHadronicInteraction: Energy ranges of two models fully overlapping");
153  }
154  rand = G4UniformRand();
155  if( emi1 < emi2 )
156  {
157  if( (ema1-kineticEnergy) < rand*(ema1-emi2) ) {
158  mem = memor2;
159  } else {
160  mem = memory;
161  }
162  } else {
163  if( (ema2-kineticEnergy) < rand*(ema2-emi1) ) {
164  mem = memory;
165  } else {
166  mem = memor2;
167  }
168  }
169  break;
170  default:
171  throw G4HadronicException(__FILE__, __LINE__,
172  "GetHadronicInteraction: More than two competing models in this energy range");
173  }
174 
175  return theHadronicInteraction[mem];
176 }
177 
178 
181  const G4Material* aMaterial,
182  const G4Element* anElement) const
183 {
185  throw G4HadronicException(__FILE__, __LINE__,
186  "GetHadronicInteraction: NO MODELS STORED");
187  }
188  G4int cou = 0, memory = 0, memor2 = 0;
189  G4double emi1 = 0.0, ema1 = 0.0, emi2 = 0.0, ema2 = 0.0;
190 
191  for (G4int i = 0; i<theHadronicInteractionCounter; ++i) {
192  G4double low = theHadronicInteraction[i]->GetMinEnergy( aMaterial, anElement );
193  // Work-around for particles with 0 kinetic energy, which still
194  // require a model to return a ParticleChange
195  //if (low == 0.) low = -DBL_MIN;
196  G4double high = theHadronicInteraction[i]->GetMaxEnergy( aMaterial, anElement );
197  if (low <= kineticEnergy && high > kineticEnergy) {
198  ++cou;
199  emi2 = emi1;
200  ema2 = ema1;
201  emi1 = low;
202  ema1 = high;
203  memor2 = memory;
204  memory = i;
205  }
206  }
207 
208  G4int mem = -1;
209  G4double rand;
210  switch (cou) {
211  case 0:
212  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="
213  <<theHadronicInteractionCounter<<", Ek="
214  <<kineticEnergy<<", Material = "<<aMaterial->GetName()
215  <<", Element = "
216  <<anElement->GetName()<<G4endl;
217  for( G4int j=0; j<theHadronicInteractionCounter; ++j)
218  {
220  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
221  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
222  }
223  throw G4HadronicException(__FILE__, __LINE__,
224  "GetHadronicInteraction: No Model found");
225  return 0;
226  case 1:
227  mem = memory;
228  break;
229  case 2:
230  if( (emi2<=emi1 && ema2>=ema1) || (emi2>=emi1 && ema2<=ema1) )
231  {
232  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="
233  <<theHadronicInteractionCounter<<", Ek="
234  <<kineticEnergy<<", Material = "<<aMaterial->GetName()
235  <<", Element = "
236  <<anElement->GetName()<<G4endl;
237  for( G4int j=0; j<theHadronicInteractionCounter; ++j)
238  {
240  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
241  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
242  }
243  throw G4HadronicException(__FILE__, __LINE__,
244  "GetHadronicInteraction: Energy ranges of two models fully overlapping");
245  }
246  rand = G4UniformRand();
247  if( emi1 < emi2 )
248  {
249  if( (ema1-kineticEnergy) < rand*(ema1-emi2) ) {
250  mem = memor2;
251  } else {
252  mem = memory;
253  }
254  } else {
255  if( (ema2-kineticEnergy) < rand*(ema2-emi1) ) {
256  mem = memory;
257  } else {
258  mem = memor2;
259  }
260  }
261  break;
262  default:
263  throw G4HadronicException(__FILE__, __LINE__,
264  "GetHadronicInteraction: More than two competing models in this energy range");
265  }
266 
267  return theHadronicInteraction[mem];
268 }
269 
270 std::vector<G4HadronicInteraction*>&
272 {
273  return theHadronicInteraction;
274 }
275 
276 #include "G4SystemOfUnits.hh"
278 {
279  G4cout << "G4EnergyRangeManager " << this << G4endl;
280  for (G4int i = 0 ; i < theHadronicInteractionCounter; i++) {
281  G4cout << " HadronicModel " << i <<":"
282  << theHadronicInteraction[i]->GetModelName() << G4endl;
283  if (verbose > 0) {
284  G4cout << " Minimum Energy "
285  << theHadronicInteraction[i]->GetMinEnergy()/GeV << " [GeV], "
286  << "Maximum Energy "
287  << theHadronicInteraction[i]->GetMaxEnergy()/GeV << " [GeV]"
288  << G4endl;
289  }
290  }
291 }
292 
293 void
295 {
296  for ( std::vector<G4HadronicInteraction*>::iterator
297  it = theHadronicInteraction.begin() ; it != theHadronicInteraction.end() ; it++ ) {
298  (*it)->BuildPhysicsTable( aParticleType );
299  }
300 }
301  /* end of file */
302 
G4double GetMinEnergy() const
void RegisterMe(G4HadronicInteraction *a)
G4double GetMaxEnergy() const
void Dump(G4int verbose=0)
G4HadronicInteraction * GetHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement) const
const G4String & GetName() const
Definition: G4Material.hh:178
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
std::vector< G4HadronicInteraction * > theHadronicInteraction
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
G4EnergyRangeManager & operator=(const G4EnergyRangeManager &right)
G4double GetKineticEnergy() const
static const double GeV
Definition: G4SIunits.hh:214
void BuildPhysicsTable(const G4ParticleDefinition &)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127