Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 71734 2013-06-21 08:53:11Z 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 
39 
41  : theHadronicInteractionCounter(0)
42 {
43  for (G4int i = 0; i < G4EnergyRangeManager::MAX_NUMBER_OF_MODELS; i++)
44  theHadronicInteraction[i] = 0;
45 }
46 
47 
49 {
50  if (this != &right) {
51  theHadronicInteractionCounter = right.theHadronicInteractionCounter;
52  for (G4int i = 0; i < theHadronicInteractionCounter; ++i)
53  theHadronicInteraction[i] = right.theHadronicInteraction[i];
54  }
55 }
56 
57 
60 {
61  if (this != &right) {
62  theHadronicInteractionCounter = right.theHadronicInteractionCounter;
63  for (G4int i=0; i<theHadronicInteractionCounter; ++i)
64  theHadronicInteraction[i] = right.theHadronicInteraction[i];
65  }
66  return *this;
67 }
68 
69 
71 {
72  if (theHadronicInteractionCounter+1 > MAX_NUMBER_OF_MODELS) {
73  throw G4HadronicException(__FILE__, __LINE__,"RegisterMe: TOO MANY MODELS");
74  }
75  theHadronicInteraction[ theHadronicInteractionCounter++ ] = a;
76 }
77 
78 
81  const G4Material* aMaterial,
82  const G4Element* anElement) const
83 {
85  if (counter == 0) throw G4HadronicException(__FILE__, __LINE__,
86  "GetHadronicInteraction: NO MODELS STORED");
87 
88  G4int cou = 0, memory = 0, memor2 = 0;
89  G4double emi1 = 0.0, ema1 = 0.0, emi2 = 0.0, ema2 = 0.0;
90 
91  for (G4int i = 0; i < counter; i++) {
92  G4double low = theHadronicInteraction[i]->GetMinEnergy( aMaterial, anElement );
93  // Work-around for particles with 0 kinetic energy, which still
94  // require a model to return a ParticleChange
95  if (low == 0.) low = -DBL_MIN;
96  G4double high = theHadronicInteraction[i]->GetMaxEnergy( aMaterial, anElement );
97  if (low < kineticEnergy && high >= kineticEnergy) {
98  ++cou;
99  emi2 = emi1;
100  ema2 = ema1;
101  emi1 = low;
102  ema1 = high;
103  memor2 = memory;
104  memory = i;
105  }
106  }
107 
108  G4int mem = -1;
109  G4double rand;
110  switch (cou) {
111  case 0:
112  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="<<counter<<", Ek="
113  <<kineticEnergy<<", Material = "<<aMaterial->GetName()<<", Element = "
114  <<anElement->GetName()<<G4endl;
115  for( G4int j=0; j<counter; j++ )
116  {
117  G4HadronicInteraction* HInt=theHadronicInteraction[j];
118  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
119  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
120  }
121  throw G4HadronicException(__FILE__, __LINE__,
122  "GetHadronicInteraction: No Model found");
123  return 0;
124  case 1:
125  mem = memory;
126  break;
127  case 2:
128  if( (emi2<=emi1 && ema2>=ema1) || (emi2>=emi1 && ema2<=ema1) )
129  {
130  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="<<counter<<", Ek="
131  <<kineticEnergy<<", Material = "<<aMaterial->GetName()<<", Element = "
132  <<anElement->GetName()<<G4endl;
133  if(counter) for( G4int j=0; j<counter; j++ )
134  {
135  G4HadronicInteraction* HInt=theHadronicInteraction[j];
136  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
137  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
138  }
139  throw G4HadronicException(__FILE__, __LINE__,
140  "GetHadronicInteraction: Energy ranges of two models fully overlapping");
141  }
142  rand = G4UniformRand();
143  if( emi1 < emi2 )
144  {
145  if( (ema1-kineticEnergy)/(ema1-emi2)<rand )
146  mem = memor2;
147  else
148  mem = memory;
149  } else {
150  if( (ema2-kineticEnergy)/(ema2-emi1)<rand )
151  mem = memory;
152  else
153  mem = memor2;
154  }
155  break;
156  default:
157  throw G4HadronicException(__FILE__, __LINE__,
158  "GetHadronicInteraction: More than two competing models in this energy range");
159  }
160 
161  return theHadronicInteraction[mem];
162 }
163 
164 
165 #include "G4SystemOfUnits.hh"
167 {
168  G4cout << "G4EnergyRangeManager " << this << G4endl;
169  for (G4int i = 0 ; i < theHadronicInteractionCounter; i++) {
170  G4cout << " HadronicModel " << i <<":"
171  << theHadronicInteraction[i]->GetModelName() << G4endl;
172  if (verbose > 0) {
173  G4cout << " Minimum Energy " << theHadronicInteraction[i]->GetMinEnergy()/GeV << " [GeV], "
174  << "Maximum Energy " << theHadronicInteraction[i]->GetMaxEnergy()/GeV << " [GeV]"
175  << G4endl;
176  }
177  }
178 }
179  /* end of file */
180 
G4double GetMinEnergy() const
void RegisterMe(G4HadronicInteraction *a)
G4double GetMaxEnergy() const
void Dump(G4int verbose=0)
const G4String & GetName() const
Definition: G4Material.hh:176
const G4String & GetModelName() const
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4EnergyRangeManager & operator=(const G4EnergyRangeManager &right)
G4int GetHadronicInteractionCounter() const
#define DBL_MIN
Definition: templates.hh:75
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
G4HadronicInteraction * GetHadronicInteraction(const G4double kineticEnergy, const G4Material *aMaterial, const G4Element *anElement) const