source: trunk/source/processes/hadronic/management/src/G4EnergyRangeManager.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

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