G4GammaConversion.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 //------------------ G4GammaConversion physics process -------------------------
00030 //                   by Michel Maire, 24 May 1996
00031 //
00032 // 11-06-96 Added SelectRandomAtom() method, M.Maire
00033 // 21-06-96 SetCuts implementation, M.Maire
00034 // 24-06-96 simplification in ComputeCrossSectionPerAtom, M.Maire
00035 // 24-06-96 in DoIt : change the particleType stuff, M.Maire
00036 // 25-06-96 modification in the generation of the teta angle, M.Maire
00037 // 16-09-96 minors optimisations in DoIt. Thanks to P.Urban
00038 //          dynamical array PartialSumSigma
00039 // 13-12-96 fast sampling of epsil below 2 MeV, L.Urban
00040 // 14-01-97 crossection table + meanfreepath table.
00041 //          PartialSumSigma removed, M.Maire
00042 // 14-01-97 in DoIt the positron is always created, even with Ekine=0,
00043 //          for further annihilation, M.Maire
00044 // 14-03-97 new Physics scheme for geant4alpha, M.Maire
00045 // 28-03-97 protection in BuildPhysicsTable, M.Maire
00046 // 19-06-97 correction in ComputeCrossSectionPerAtom, L.Urban
00047 // 04-06-98 in DoIt, secondary production condition:
00048 //            range>std::min(threshold,safety)
00049 // 13-08-98 new methods SetBining() PrintInfo()
00050 // 28-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
00051 // 11-07-01 PostStepDoIt - sampling epsil: power(rndm,0.333333)
00052 // 13-07-01 DoIt: suppression of production cut for the (e-,e+) (mma)
00053 // 06-08-01 new methods Store/Retrieve PhysicsTable (mma)
00054 // 06-08-01 BuildThePhysicsTable() called from constructor (mma)
00055 // 17-09-01 migration of Materials to pure STL (mma)
00056 // 20-09-01 DoIt: fminimalEnergy = 1*eV (mma)
00057 // 01-10-01 come back to BuildPhysicsTable(const G4ParticleDefinition&)
00058 // 11-01-02 ComputeCrossSection: correction of extrapolation below EnergyLimit
00059 // 21-03-02 DoIt: correction of the e+e- angular distribution (bug 363) mma
00060 // 08-11-04 Remove of Store/Retrieve tables (V.Ivantchenko)
00061 // 19-04-05 Migrate to model interface and inherit 
00062 //          from G4VEmProcess (V.Ivanchenko) 
00063 // 04-05-05, Make class to be default (V.Ivanchenko)
00064 // 09-08-06, add SetModel(G4VEmModel*) (mma)
00065 // 12-09-06, move SetModel(G4VEmModel*) in G4VEmProcess (mma)
00066 // -----------------------------------------------------------------------------
00067 
00068 #include "G4GammaConversion.hh"
00069 #include "G4PhysicalConstants.hh"
00070 #include "G4SystemOfUnits.hh"
00071 #include "G4BetheHeitlerModel.hh"
00072 #include "G4PairProductionRelModel.hh"
00073 #include "G4Electron.hh"
00074 
00075 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00076 
00077 using namespace std;
00078 
00079 G4GammaConversion::G4GammaConversion(const G4String& processName,
00080   G4ProcessType type):G4VEmProcess (processName, type),
00081     isInitialised(false)
00082 {
00083   SetMinKinEnergy(2.0*electron_mass_c2);
00084   SetProcessSubType(fGammaConversion);
00085   SetStartFromNullFlag(true);
00086   SetBuildTableFlag(true);
00087   SetSecondaryParticle(G4Electron::Electron());
00088 }
00089 
00090 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00091  
00092 G4GammaConversion::~G4GammaConversion()
00093 {}
00094 
00095 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00096 
00097 G4bool G4GammaConversion::IsApplicable(const G4ParticleDefinition& p)
00098 {
00099   return (&p == G4Gamma::Gamma());
00100 }
00101 
00102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00103 
00104 void G4GammaConversion::InitialiseProcess(const G4ParticleDefinition*)
00105 {
00106   if(!isInitialised) {
00107     isInitialised = true;
00108     const G4double limit = 80*GeV;
00109     G4double emin = std::max(MinKinEnergy(), 2*electron_mass_c2);
00110     G4double emax = MaxKinEnergy();
00111     SetMinKinEnergy(emin);
00112 
00113     if(!EmModel(1)) { SetEmModel(new G4BetheHeitlerModel(), 1); }
00114     EmModel(1)->SetLowEnergyLimit(emin);
00115     G4double ehigh = std::min(emax,limit);
00116     ehigh = std::min(ehigh,EmModel(1)->HighEnergyLimit());
00117     EmModel(1)->SetHighEnergyLimit(ehigh);
00118     AddEmModel(1, EmModel(1));
00119 
00120     if(emax > ehigh) {
00121       if(!EmModel(2)) { SetEmModel(new G4PairProductionRelModel(), 2); }
00122       EmModel(2)->SetLowEnergyLimit(ehigh);
00123       EmModel(2)->SetHighEnergyLimit(emax);
00124       AddEmModel(2, EmModel(2));
00125     }
00126   } 
00127 }
00128 
00129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00130 
00131 G4double G4GammaConversion::MinPrimaryEnergy(const G4ParticleDefinition*,
00132                                              const G4Material*)
00133 {
00134   return 2*electron_mass_c2;
00135 }
00136 
00137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00138 
00139 void G4GammaConversion::PrintInfo()
00140 {}         
00141 
00142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

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