Geant4  9.6.p02
 All Classes 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 //
27 // $Id$
28 //
29  // Hadronic Process: Energy Range Manager
30  // original by H.P. Wellisch
31  // modified by J.L. Chuma, TRIUMF, 22-Nov-1996
32  // Last modified: 24-Mar-1997
33  // fix in the counter-hndling: H.P. Wellisch 04-Apr-97
34  // throw an exception if no model found: J.L. Chuma 04-Apr-97
35 
36 #include "G4EnergyRangeManager.hh"
37 #include "Randomize.hh"
38 #include "G4HadronicException.hh"
39 
40 
42  : theHadronicInteractionCounter(0)
43 {
44  for (G4int i = 0; i < G4EnergyRangeManager::MAX_NUMBER_OF_MODELS; i++)
45  theHadronicInteraction[i] = 0;
46 }
47 
48 
50 {
51  if (this != &right) {
52  for (G4int i=0; i<theHadronicInteractionCounter; ++i)
53  theHadronicInteraction[i] = right.theHadronicInteraction[i];
54  theHadronicInteractionCounter = right.theHadronicInteractionCounter;
55  }
56 }
57 
58 
61 {
62  if (this != &right) {
63  for (G4int i=0; i<theHadronicInteractionCounter; ++i)
64  theHadronicInteraction[i] = right.theHadronicInteraction[i];
65  theHadronicInteractionCounter = right.theHadronicInteractionCounter;
66  }
67  return *this;
68 }
69 
70 
72 {
73  if( theHadronicInteractionCounter+1 > MAX_NUMBER_OF_MODELS )
74  {
75  throw G4HadronicException(__FILE__, __LINE__,"RegisterMe: TOO MANY MODELS");
76  }
77  theHadronicInteraction[ theHadronicInteractionCounter++ ] = a;
78 }
79 
80 
83  const G4double kineticEnergy,
84  const G4Material *aMaterial,
85  const G4Element *anElement ) const
86  {
88  if( counter == 0 )
89  throw G4HadronicException(__FILE__, __LINE__,
90  "GetHadronicInteraction: NO MODELS STORED");
91 
92  G4int cou = 0, memory = 0, memor2 = 0;
93  G4double emi1 = 0.0, ema1 = 0.0, emi2 = 0.0, ema2 = 0.0;
94  for( G4int i=0; i<counter; i++ )
95  {
96  G4double low = theHadronicInteraction[i]->GetMinEnergy( aMaterial, anElement );
97  // Work-around for particles with 0 kinetic energy, which still
98  // require a model to return a ParticleChange
99  if (low == 0.) low = -DBL_MIN;
100  G4double high = theHadronicInteraction[i]->GetMaxEnergy( aMaterial, anElement );
101  if( low < kineticEnergy && high >= kineticEnergy )
102  {
103  ++cou;
104  emi2 = emi1;
105  ema2 = ema1;
106  emi1 = low;
107  ema1 = high;
108  memor2 = memory;
109  memory = i;
110  }
111  }
112  G4int mem=-1;
113  G4double rand;
114  switch ( cou )
115  {
116  case 0:
117  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="<<counter<<", Ek="
118  <<kineticEnergy<<", Material = "<<aMaterial->GetName()<<", Element = "
119  <<anElement->GetName()<<G4endl;
120  for( G4int j=0; j<counter; j++ )
121  {
122  G4HadronicInteraction* HInt=theHadronicInteraction[j];
123  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
124  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
125  }
126  throw G4HadronicException(__FILE__, __LINE__,
127  "GetHadronicInteraction: No Model found");
128  return 0;
129  case 1:
130  mem = memory;
131  break;
132  case 2:
133  if( (emi2<=emi1 && ema2>=ema1) || (emi2>=emi1 && ema2<=ema1) )
134  {
135  G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="<<counter<<", Ek="
136  <<kineticEnergy<<", Material = "<<aMaterial->GetName()<<", Element = "
137  <<anElement->GetName()<<G4endl;
138  if(counter) for( G4int j=0; j<counter; j++ )
139  {
140  G4HadronicInteraction* HInt=theHadronicInteraction[j];
141  G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
142  <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
143  }
144  throw G4HadronicException(__FILE__, __LINE__,
145  "GetHadronicInteraction: Energy ranges of two models fully overlapping");
146  }
147  rand = G4UniformRand();
148  if( emi1 < emi2 )
149  {
150  if( (ema1-kineticEnergy)/(ema1-emi2)<rand )
151  mem = memor2;
152  else
153  mem = memory;
154  } else {
155  if( (ema2-kineticEnergy)/(ema2-emi1)<rand )
156  mem = memory;
157  else
158  mem = memor2;
159  }
160  break;
161  default:
162  throw G4HadronicException(__FILE__, __LINE__,
163  "GetHadronicInteraction: More than two competing models in this energy range");
164  }
165  return theHadronicInteraction[mem];
166  }
167 
168  /* end of file */
169