G4EnergyRangeManager.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029  // Hadronic Process: Energy Range Manager
00030  // original by H.P. Wellisch
00031  // modified by J.L. Chuma, TRIUMF, 22-Nov-1996
00032  // Last modified: 24-Mar-1997
00033  // fix in the counter-hndling: H.P. Wellisch 04-Apr-97
00034  // throw an exception if no model found:  J.L. Chuma  04-Apr-97
00035  
00036 #include "G4EnergyRangeManager.hh"
00037 #include "Randomize.hh"
00038 #include "G4HadronicException.hh"
00039 
00040 
00041 G4EnergyRangeManager::G4EnergyRangeManager()
00042  : theHadronicInteractionCounter(0)
00043 {
00044   for (G4int i = 0; i < G4EnergyRangeManager::MAX_NUMBER_OF_MODELS; i++)
00045     theHadronicInteraction[i] = 0;
00046 }
00047 
00048 
00049 G4EnergyRangeManager::G4EnergyRangeManager(const G4EnergyRangeManager& right)
00050 {
00051   if (this != &right) {
00052     for (G4int i=0; i<theHadronicInteractionCounter; ++i)
00053       theHadronicInteraction[i] = right.theHadronicInteraction[i];
00054     theHadronicInteractionCounter = right.theHadronicInteractionCounter;
00055   }
00056 }
00057  
00058 
00059 G4EnergyRangeManager&G4EnergyRangeManager::operator=(
00060    const G4EnergyRangeManager &right )
00061 {
00062   if (this != &right) {
00063     for (G4int i=0; i<theHadronicInteractionCounter; ++i)
00064       theHadronicInteraction[i] = right.theHadronicInteraction[i];
00065     theHadronicInteractionCounter = right.theHadronicInteractionCounter;
00066   }
00067   return *this;
00068 }
00069 
00070 
00071 void G4EnergyRangeManager::RegisterMe(G4HadronicInteraction* a)
00072 {
00073     if( theHadronicInteractionCounter+1 > MAX_NUMBER_OF_MODELS )
00074     {
00075       throw G4HadronicException(__FILE__, __LINE__,"RegisterMe: TOO MANY MODELS");
00076     }
00077     theHadronicInteraction[ theHadronicInteractionCounter++ ] = a;
00078 }
00079 
00080  
00081 G4HadronicInteraction*
00082   G4EnergyRangeManager::GetHadronicInteraction(
00083    const G4double kineticEnergy,
00084    const G4Material *aMaterial,
00085    const G4Element *anElement ) const
00086   {
00087     G4int counter = GetHadronicInteractionCounter();
00088     if( counter == 0 )
00089       throw G4HadronicException(__FILE__, __LINE__,
00090                                "GetHadronicInteraction: NO MODELS STORED");
00091 
00092     G4int cou = 0, memory = 0, memor2 = 0;
00093     G4double emi1 = 0.0, ema1 = 0.0, emi2 = 0.0, ema2 = 0.0;
00094     for( G4int i=0; i<counter; i++ )
00095     {
00096       G4double low  = theHadronicInteraction[i]->GetMinEnergy( aMaterial, anElement );
00097       // Work-around for particles with 0 kinetic energy, which still
00098       // require a model to return a ParticleChange
00099       if (low == 0.) low = -DBL_MIN;
00100       G4double high = theHadronicInteraction[i]->GetMaxEnergy( aMaterial, anElement );
00101       if( low < kineticEnergy && high >= kineticEnergy )
00102       {
00103         ++cou;
00104         emi2 = emi1;
00105         ema2 = ema1;
00106         emi1 = low;
00107         ema1 = high;
00108         memor2 = memory;
00109         memory = i;
00110       }
00111     }
00112     G4int mem=-1;
00113     G4double rand;
00114     switch ( cou )
00115     {
00116      case 0:
00117        G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="<<counter<<", Ek="
00118              <<kineticEnergy<<", Material = "<<aMaterial->GetName()<<", Element = "
00119              <<anElement->GetName()<<G4endl;
00120        for( G4int j=0; j<counter; j++ )
00121        {
00122          G4HadronicInteraction* HInt=theHadronicInteraction[j];
00123          G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
00124                <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
00125        }
00126        throw G4HadronicException(__FILE__, __LINE__,
00127           "GetHadronicInteraction: No Model found");
00128        return 0;
00129      case 1:
00130        mem = memory;
00131        break;
00132      case 2:
00133        if( (emi2<=emi1 && ema2>=ema1) || (emi2>=emi1 && ema2<=ema1) )
00134        {
00135          G4cout<<"G4EnergyRangeManager:GetHadronicInteraction: counter="<<counter<<", Ek="
00136                <<kineticEnergy<<", Material = "<<aMaterial->GetName()<<", Element = "
00137                <<anElement->GetName()<<G4endl;
00138          if(counter) for( G4int j=0; j<counter; j++ )
00139          {
00140            G4HadronicInteraction* HInt=theHadronicInteraction[j];
00141            G4cout<<"*"<<j<<"* low=" <<HInt->GetMinEnergy(aMaterial,anElement)
00142                <<", high="<<HInt->GetMaxEnergy(aMaterial,anElement)<<G4endl;
00143          }
00144          throw G4HadronicException(__FILE__, __LINE__,
00145                "GetHadronicInteraction: Energy ranges of two models fully overlapping");
00146        }
00147        rand = G4UniformRand();
00148        if( emi1 < emi2 )
00149        {
00150          if( (ema1-kineticEnergy)/(ema1-emi2)<rand )
00151            mem = memor2;
00152          else
00153            mem = memory;
00154        } else {
00155          if( (ema2-kineticEnergy)/(ema2-emi1)<rand )
00156            mem = memory;
00157          else
00158            mem = memor2;
00159        }
00160        break;
00161      default:
00162       throw G4HadronicException(__FILE__, __LINE__,
00163         "GetHadronicInteraction: More than two competing models in this energy range");
00164     }
00165     return theHadronicInteraction[mem];
00166   } 
00167  
00168  /* end of file */
00169  

Generated on Mon May 27 17:48:10 2013 for Geant4 by  doxygen 1.4.7