G4PEEffectFluoModel.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 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4PEEffectFluoModel
00034 //
00035 // Author:        Vladimir Ivanchenko on base of G4PEEffectModel
00036 //
00037 // Creation date: 13.06.2010
00038 //
00039 // Modifications:
00040 //
00041 // Class Description:
00042 // Implementation of the photo-electric effect with deexcitation
00043 //
00044 // -------------------------------------------------------------------
00045 //
00046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00047 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00048 
00049 #include "G4PEEffectFluoModel.hh"
00050 #include "G4PhysicalConstants.hh"
00051 #include "G4SystemOfUnits.hh"
00052 #include "G4Electron.hh"
00053 #include "G4Gamma.hh"
00054 #include "Randomize.hh"
00055 #include "G4DataVector.hh"
00056 #include "G4ParticleChangeForGamma.hh"
00057 #include "G4VAtomDeexcitation.hh"
00058 #include "G4LossTableManager.hh"
00059 #include "G4SauterGavrilaAngularDistribution.hh"
00060 
00061 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00062 
00063 using namespace std;
00064 
00065 G4PEEffectFluoModel::G4PEEffectFluoModel(const G4String& nam)
00066   : G4VEmModel(nam)
00067 {
00068   theGamma    = G4Gamma::Gamma();
00069   theElectron = G4Electron::Electron();
00070   fminimalEnergy = 1.0*eV;
00071   SetDeexcitationFlag(true);
00072   fParticleChange = 0;
00073   fAtomDeexcitation = 0;
00074 
00075   // default generator
00076   SetAngularDistribution(new G4SauterGavrilaAngularDistribution());
00077 }
00078 
00079 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00080 
00081 G4PEEffectFluoModel::~G4PEEffectFluoModel()
00082 {}
00083 
00084 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00085 
00086 void G4PEEffectFluoModel::Initialise(const G4ParticleDefinition*,
00087                                      const G4DataVector&)
00088 {
00089   fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
00090   if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
00091 }
00092 
00093 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00094 
00095 G4double 
00096 G4PEEffectFluoModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
00097                                                 G4double energy,
00098                                                 G4double Z, G4double,
00099                                                 G4double, G4double)
00100 {
00101   G4double* SandiaCof = G4SandiaTable::GetSandiaCofPerAtom((G4int)Z, energy);
00102 
00103   G4double energy2 = energy*energy;
00104   G4double energy3 = energy*energy2;
00105   G4double energy4 = energy2*energy2;
00106 
00107   return SandiaCof[0]/energy  + SandiaCof[1]/energy2 +
00108     SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
00109 }
00110 
00111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00112 
00113 G4double 
00114 G4PEEffectFluoModel::CrossSectionPerVolume(const G4Material* material,
00115                                            const G4ParticleDefinition*,
00116                                            G4double energy,
00117                                            G4double, G4double)
00118 {
00119   G4double* SandiaCof = 
00120     material->GetSandiaTable()->GetSandiaCofForMaterial(energy);
00121                                 
00122   G4double energy2 = energy*energy;
00123   G4double energy3 = energy*energy2;
00124   G4double energy4 = energy2*energy2;
00125           
00126   return SandiaCof[0]/energy  + SandiaCof[1]/energy2 +
00127     SandiaCof[2]/energy3 + SandiaCof[3]/energy4; 
00128 }
00129 
00130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00131 
00132 void 
00133 G4PEEffectFluoModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
00134                                        const G4MaterialCutsCouple* couple,
00135                                        const G4DynamicParticle* aDynamicPhoton,
00136                                        G4double,
00137                                        G4double)
00138 {
00139   const G4Material* aMaterial = couple->GetMaterial();
00140 
00141   G4double energy = aDynamicPhoton->GetKineticEnergy();
00142 
00143   // select randomly one element constituing the material.
00144   const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
00145   
00146   //
00147   // Photo electron
00148   //
00149 
00150   // Select atomic shell
00151   G4int nShells = anElement->GetNbOfAtomicShells();
00152   G4int i = 0;  
00153   for(; i<nShells; ++i) {
00154     /*
00155     G4cout << "i= " << i << " E(eV)= " << energy/eV 
00156            << " Eb(eV)= " << anElement->GetAtomicShell(i)/eV
00157            << "  " << anElement->GetName() 
00158            << G4endl;
00159     */
00160     if(energy >= anElement->GetAtomicShell(i)) { break; }
00161   }
00162 
00163   G4double edep = energy;
00164 
00165   // Normally one shell is available 
00166   if (i < nShells) { 
00167   
00168     G4double bindingEnergy  = anElement->GetAtomicShell(i);
00169     G4double elecKineEnergy = energy - bindingEnergy;
00170 
00171     // create photo electron
00172     //
00173     if (elecKineEnergy > fminimalEnergy) {
00174       edep = bindingEnergy;
00175       G4ThreeVector elecDirection =
00176         GetAngularDistribution()->SampleDirection(aDynamicPhoton, 
00177                                                   elecKineEnergy,
00178                                                   i, 
00179                                                   couple->GetMaterial());
00180    
00181       G4DynamicParticle* aParticle = 
00182         new G4DynamicParticle(theElectron, elecDirection, elecKineEnergy);
00183       fvect->push_back(aParticle);
00184     } 
00185 
00186     // sample deexcitation
00187     //
00188     if(fAtomDeexcitation) {
00189       G4int index = couple->GetIndex();
00190       if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
00191         G4int Z = G4lrint(anElement->GetZ());
00192         G4AtomicShellEnumerator as = G4AtomicShellEnumerator(i);
00193         const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
00194         size_t nbefore = fvect->size();
00195         fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
00196         size_t nafter = fvect->size();
00197         if(nafter > nbefore) {
00198           for (size_t j=nbefore; j<nafter; ++j) {
00199             edep -= ((*fvect)[j])->GetKineticEnergy();
00200           } 
00201         }
00202       }
00203     }
00204   }
00205 
00206   // kill primary photon
00207   fParticleChange->SetProposedKineticEnergy(0.);
00208   fParticleChange->ProposeTrackStatus(fStopAndKill);
00209   if(edep > 0.0) {
00210     fParticleChange->ProposeLocalEnergyDeposit(edep);
00211   }
00212 }
00213 
00214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

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