G4PenelopeIonisationCrossSection.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 // $Id$
00027 //
00028 // Author: Luciano Pandola
00029 //
00030 // History:
00031 // --------
00032 // 14 Mar 2012   L Pandola    First complete implementation for e-
00033 //
00034 //
00037 //
00038 #include "G4PenelopeOscillatorManager.hh"
00039 #include "G4SystemOfUnits.hh"
00040 #include "G4PenelopeIonisationCrossSection.hh"
00041 #include "G4PenelopeIonisationXSHandler.hh"
00042 #include "G4PenelopeCrossSection.hh"
00043 #include "G4Electron.hh"
00044 #include "G4AtomicTransitionManager.hh"
00045 
00046 G4PenelopeIonisationCrossSection::G4PenelopeIonisationCrossSection() : 
00047   G4VhShellCrossSection("Penelope"),shellIDTable(0),
00048   theCrossSectionHandler(0)
00049 { 
00050   oscManager = G4PenelopeOscillatorManager::GetOscillatorManager();
00051   nMaxLevels = 9;
00052 
00053   // Verbosity scale:
00054   // 0 = nothing
00055   // 1 = calculation of cross sections, file openings, sampling of atoms
00056   // 2 = entering in methods
00057   verboseLevel = 0;
00058 
00059   fLowEnergyLimit  = 10.0*eV;
00060   fHighEnergyLimit = 100.0*GeV;
00061 
00062   transitionManager = G4AtomicTransitionManager::Instance();
00063 }
00064 
00065 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00066  
00067 G4PenelopeIonisationCrossSection::~G4PenelopeIonisationCrossSection()
00068 {
00069   if (theCrossSectionHandler)
00070     delete theCrossSectionHandler;
00071 }
00072 
00073 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00074 
00075 G4double G4PenelopeIonisationCrossSection::CrossSection(G4int Z,
00076                                                         G4AtomicShellEnumerator shell,
00077                                                         G4double incidentEnergy,
00078                                                         G4double ,
00079                                                         const G4Material* material)
00080 {
00081   if (verboseLevel > 1)    
00082     G4cout << "Entering in method G4PenelopeIonisationCrossSection::CrossSection()" << G4endl;
00083         
00084   G4double cross = 0.;
00085 
00086   //Material pointer is not available
00087   if (!material)
00088     {
00089       //CRASH!
00090       G4ExceptionDescription ed;
00091       ed << "The method has been called with a NULL G4Material pointer" << G4endl;
00092       G4Exception("G4PenelopeIonisationCrossSection::CrossSection()","em2042",
00093                   FatalException,ed);
00094       return cross;
00095     }
00096 
00097   if (!theCrossSectionHandler)
00098     theCrossSectionHandler = new G4PenelopeIonisationXSHandler();
00099 
00100   G4int nmax = std::min(nMaxLevels,transitionManager->NumberOfShells(Z));
00101 
00102   if(G4int(shell) < nmax && 
00103      incidentEnergy >= fLowEnergyLimit && incidentEnergy <= fHighEnergyLimit) 
00104     {
00105       //The shells in Penelope are organized per *material*, rather than per 
00106       //element, so given a material one has to find the proper index for the 
00107       //given Z and shellID. An appropriate lookup table is used to avoid 
00108       //recalculation.
00109       G4int index = FindShellIDIndex(material,Z,shell);
00110 
00111       //Index is not available!
00112       if (index < 0)
00113         return cross;
00114 
00115       G4PenelopeCrossSection* theXS = 
00116         theCrossSectionHandler->GetCrossSectionTableForCouple(G4Electron::Electron(),
00117                                                               material,
00118                                                               0.);
00119 
00120       //Cross check that everything is fine:
00121       G4PenelopeOscillator* theOsc = oscManager->GetOscillatorIonisation(material,index);
00122       if (theOsc->GetParentZ() != Z || theOsc->GetShellFlag()-1 != G4int(shell))
00123         {
00124           //something went wrong!
00125           G4ExceptionDescription ed;
00126           ed << "There is something wrong here: it looks like the index is wrong" << G4endl;
00127           ed << "Requested: shell " << G4int(shell) << " and Z = " << Z << G4endl;
00128           ed << "Retrieved: " << theOsc->GetShellFlag()-1 << " and Z = " << theOsc->GetParentZ() << G4endl;       
00129           G4Exception("G4PenelopeIonisationCrossSection::CrossSection()","em2043",
00130                       JustWarning,ed);
00131           return cross;
00132         }
00133 
00134 
00135       G4double crossPerMolecule = (theXS) ? theXS->GetShellCrossSection(index,incidentEnergy) : 0.;
00136 
00137       //Now it must be converted in cross section per atom. I need the number of 
00138       //atoms of the given Z per molecule.
00139       G4double atomsPerMolec = oscManager->GetNumberOfZAtomsPerMolecule(material,Z);
00140       if (atomsPerMolec)
00141         cross = crossPerMolecule/atomsPerMolec;
00142      
00143 
00144       if (verboseLevel > 0)
00145         {
00146           G4cout << "Cross section of shell " << G4int(shell) << " and Z= " << Z;
00147           G4cout << " of material: " << material->GetName() << " and energy = " << incidentEnergy/keV << " keV" << G4endl;
00148           G4cout << "--> " << cross/barn << " barn" << G4endl;
00149           G4cout << "Shell binding energy: " << theOsc->GetIonisationEnergy()/eV << " eV;" ;
00150           G4cout << " resonance energy: " << theOsc->GetResonanceEnergy()/eV << "eV" << G4endl;
00151            if (verboseLevel > 2)
00152              {
00153                G4cout << "Cross section per molecule: " << crossPerMolecule/barn << " barn" << G4endl;
00154                G4cout << "Atoms " << Z << " per molecule: " << atomsPerMolec << G4endl;
00155              }
00156         }     
00157     }  
00158   
00159   return cross;
00160 }
00161 
00162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00163 std::vector<G4double> 
00164 G4PenelopeIonisationCrossSection::GetCrossSection(G4int Z,
00165                                                   G4double kinEnergy,
00166                                                   G4double, G4double, 
00167                                                   const G4Material* mat)
00168 {
00169   G4int nmax = std::min(nMaxLevels,transitionManager->NumberOfShells(Z));
00170   std::vector<G4double> vec(nmax,0.0); 
00171   for(G4int i=0; i<nmax; ++i) {
00172     vec[i] = CrossSection(Z, G4AtomicShellEnumerator(i), kinEnergy,0.,mat); 
00173   }
00174   return vec;
00175 }
00176 
00177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00178 
00179 std::vector<G4double> 
00180 G4PenelopeIonisationCrossSection::Probabilities(G4int Z,
00181                                                  G4double kinEnergy,
00182                                                  G4double,
00183                                                  G4double,
00184                                                  const G4Material* mat)
00185 {
00186   std::vector<G4double> vec = GetCrossSection(Z, kinEnergy,0,0,mat);
00187   size_t n = vec.size();
00188   size_t i=0;
00189   G4double sum = 0.0;
00190   for(i=0; i<n; ++i) { sum += vec[i]; }
00191   if(sum > 0.0) { 
00192     sum = 1.0/sum; 
00193     for(i=0; i<n; ++i) { vec[i] = vec[i]*sum; }
00194   }
00195   return vec;
00196 }
00197 
00198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00199 
00200 G4int G4PenelopeIonisationCrossSection::FindShellIDIndex(const G4Material* mat,
00201                                                          G4int Z,
00202                                                          G4AtomicShellEnumerator shell)
00203 {
00204    if (verboseLevel > 1)    
00205      G4cout << "Entering in method G4PenelopeIonisationCrossSection::FindShellIDIndex()" << G4endl;
00206 
00207   if (!shellIDTable)
00208     shellIDTable = new std::map< std::pair<const G4Material*,G4int>, G4DataVector*>;
00209  
00210   std::pair<const G4Material*,G4int> theKey = std::make_pair(mat,Z);
00211   G4int result = -1;
00212   G4int ishell = G4int(shell);
00213   
00214   if (shellIDTable->count(theKey)) //table already built, and containing the element
00215     {   
00216       if (verboseLevel > 2)
00217         G4cout << "FindShellIDIndex: Table already built for " << mat->GetName() << G4endl;
00218       G4DataVector* theVec = shellIDTable->find(theKey)->second;
00219          
00220       if (ishell>=0 && ishell < (G4int) theVec->size()) //check we are not off-boundary
00221         result = (G4int) (*theVec)[ishell];    
00222       else
00223         {
00224           G4ExceptionDescription ed;
00225           ed << "Shell ID: " << ishell << " not available for material " << mat->GetName() << " and Z = " << 
00226             Z << G4endl;
00227           G4Exception("G4PenelopeIonisationCrossSection::FindShellIDIndex()","em2041",JustWarning,
00228                       ed);
00229           return -1;
00230         }          
00231     }
00232   else
00233     {
00234       if (verboseLevel > 2)
00235         G4cout << "FindShellIDIndex: Table to be built for " << mat->GetName() << G4endl;
00236       //Not contained: look for it
00237       G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableIonisation(mat);
00238       size_t numberOfOscillators = theTable->size();
00239           
00240       //oscillator loop
00241       //initialize everything at -1
00242       G4DataVector* dat = new G4DataVector(nMaxLevels,-1);
00243       for (size_t iosc=0;iosc<numberOfOscillators;iosc++)
00244         {      
00245           G4PenelopeOscillator* theOsc = (*theTable)[iosc];
00246           //level is found!
00247           if (theOsc->GetParentZ() == Z)
00248             {
00249               //individual shells relative to the given material            
00250               G4int shFlag = theOsc->GetShellFlag();      
00251               //Notice: GetShellFlag() starts from 1, the G4AtomicShellEnumerator from 0
00252               if (shFlag < 30)
00253                 (*dat)[shFlag-1] = (G4double) iosc; //index of the given shell
00254               if ((shFlag-1) == ishell) // this is what we were looking for
00255                 result = (G4int) iosc;
00256             }      
00257         }
00258       shellIDTable->insert(std::make_pair(theKey,dat));
00259     }
00260 
00261   if (verboseLevel > 1)    
00262     G4cout << "Leaving method G4PenelopeIonisationCrossSection::FindShellIDIndex() with index = " << result << G4endl;
00263 
00264   return result;
00265   
00266 }

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