// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // $Id: G4MuPairProductionModel.cc,v 1.39 2008/07/22 16:11:34 vnivanch Exp $ // GEANT4 tag $Name: geant4-09-02 $ // // ------------------------------------------------------------------- // // GEANT4 Class file // // // File name: G4MuPairProductionModel // // Author: Vladimir Ivanchenko on base of Laszlo Urban code // // Creation date: 24.06.2002 // // Modifications: // // 04-12-02 Change G4DynamicParticle constructor in PostStep (V.Ivanchenko) // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko) // 24-01-03 Fix for compounds (V.Ivanchenko) // 27-01-03 Make models region aware (V.Ivanchenko) // 13-02-03 Add model (V.Ivanchenko) // 06-06-03 Fix in cross section calculation for high energy (V.Ivanchenko) // 20-10-03 2*xi in ComputeDDMicroscopicCrossSection (R.Kokoulin) // 8 integration points in ComputeDMicroscopicCrossSection // 12-01-04 Take min cut of e- and e+ not its sum (V.Ivanchenko) // 10-02-04 Update parameterisation using R.Kokoulin model (V.Ivanchenko) // 28-04-04 For complex materials repeat calculation of max energy for each // material (V.Ivanchenko) // 01-11-04 Fix bug inside ComputeDMicroscopicCrossSection (R.Kokoulin) // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) // 03-08-05 Add SetParticle method (V.Ivantchenko) // 23-10-05 Add protection in sampling of e+e- pair energy needed for // low cuts (V.Ivantchenko) // 13-02-06 Add ComputeCrossSectionPerAtom (mma) // 24-04-07 Add protection in SelectRandomAtom method (V.Ivantchenko) // 12-05-06 Updated sampling (use cut) in SelectRandomAtom (A.Bogdanov) // 11-10-07 Add ignoreCut flag (V.Ivanchenko) // // Class Description: // // // ------------------------------------------------------------------- // //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... #include "G4MuPairProductionModel.hh" #include "G4Electron.hh" #include "G4Positron.hh" #include "G4MuonMinus.hh" #include "G4MuonPlus.hh" #include "Randomize.hh" #include "G4Material.hh" #include "G4Element.hh" #include "G4ElementVector.hh" #include "G4ProductionCutsTable.hh" #include "G4ParticleChangeForLoss.hh" #include "G4ParticleChangeForGamma.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... // static members // G4double G4MuPairProductionModel::zdat[]={1., 4., 13., 29., 92.}; G4double G4MuPairProductionModel::adat[]={1.01, 9.01, 26.98, 63.55, 238.03}; G4double G4MuPairProductionModel::tdat[]={1.e3, 1.e4, 1.e5, 1.e6, 1.e7, 1.e8, 1.e9, 1.e10}; G4double G4MuPairProductionModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917, 0.7628, 0.8983, 0.9801 }; G4double G4MuPairProductionModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813, 0.1569, 0.1112, 0.0506 }; //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... using namespace std; G4MuPairProductionModel::G4MuPairProductionModel(const G4ParticleDefinition* p, const G4String& nam) : G4VEmModel(nam), particle(0), factorForCross(4.*fine_structure_const*fine_structure_const *classic_electr_radius*classic_electr_radius/(3.*pi)), sqrte(sqrt(exp(1.))), currentZ(0), fParticleChange(0), minPairEnergy(4.*electron_mass_c2), lowestKinEnergy(1.*GeV), nzdat(5), ntdat(8), nbiny(1000), nmaxElements(0), ymin(-5.), ymax(0.), dy((ymax-ymin)/nbiny), samplingTablesAreFilled(false) { SetLowEnergyLimit(minPairEnergy); nist = G4NistManager::Instance(); theElectron = G4Electron::Electron(); thePositron = G4Positron::Positron(); if(p) SetParticle(p); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4MuPairProductionModel::~G4MuPairProductionModel() {} //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void G4MuPairProductionModel::Initialise(const G4ParticleDefinition* p, const G4DataVector&) { if (!samplingTablesAreFilled) { if(p) SetParticle(p); MakeSamplingTables(); } if(!fParticleChange) { if(pParticleChange) fParticleChange = reinterpret_cast(pParticleChange); else fParticleChange = new G4ParticleChangeForLoss(); } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4double G4MuPairProductionModel::ComputeDEDXPerVolume( const G4Material* material, const G4ParticleDefinition*, G4double kineticEnergy, G4double cutEnergy) { G4double dedx = 0.0; if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy) return dedx; const G4ElementVector* theElementVector = material->GetElementVector(); const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector(); // loop for elements in the material for (size_t i=0; iGetNumberOfElements(); i++) { G4double Z = (*theElementVector)[i]->GetZ(); SetCurrentElement(Z); G4double tmax = MaxSecondaryEnergy(particle, kineticEnergy); G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax); dedx += loss*theAtomicNumDensityVector[i]; } if (dedx < 0.) dedx = 0.; return dedx; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4double G4MuPairProductionModel::ComputMuPairLoss(G4double Z, G4double tkin, G4double cutEnergy, G4double tmax) { SetCurrentElement(Z); G4double loss = 0.0; G4double cut = std::min(cutEnergy,tmax); if(cut <= minPairEnergy) return loss; // calculate the rectricted loss // numerical integration in log(PairEnergy) G4double ak1=6.9; G4double ak2=1.0; G4double aaa = log(minPairEnergy); G4double bbb = log(cut); G4int kkk = (G4int)((bbb-aaa)/ak1+ak2); if (kkk > 8) kkk = 8; G4double hhh = (bbb-aaa)/(G4double)kkk; G4double x = aaa; for (G4int l=0 ; l 8) kkk = 8; G4double hhh = (bbb-aaa)/float(kkk); G4double x = aaa; for(G4int l=0; l 0.) { G4double zeta2 = 0.058*log(totalEnergy/(particleMass+g2*z13*totalEnergy))-0.14; zeta = zeta1/zeta2 ; } G4double z2 = Z*(Z+zeta); G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy); G4double a0 = totalEnergy*residEnergy; G4double a1 = pairEnergy*pairEnergy/a0; G4double bet = 0.5*a1; G4double xi0 = 0.25*massratio2*a1; G4double del = c8/a0; G4double rta3 = sqrt(a3); G4double tmnexp = alf/(1. + rta3) + del*rta3; if(tmnexp >= 1.0) return cross; G4double tmn = log(tmnexp); G4double sum = 0.; // Gaussian integration in ln(1-ro) ( with 8 points) for (G4int i=0; i<8; i++) { G4double a4 = exp(tmn*xgi[i]); // a4 = (1.-asymmetry) G4double a5 = a4*(2.-a4) ; G4double a6 = 1.-a5 ; G4double a7 = 1.+a6 ; G4double a9 = 3.+a6 ; G4double xi = xi0*a5 ; G4double xii = 1./xi ; G4double xi1 = 1.+xi ; G4double screen = screen0*xi1/a5 ; G4double yeu = 5.-a6+4.*bet*a7 ; G4double yed = 2.*(1.+3.*bet)*log(3.+xii)-a6-a1*(2.-a6) ; G4double ye1 = 1.+yeu/yed ; G4double ale=log(bbb/z13*sqrt(xi1*ye1)/(1.+screen*ye1)) ; G4double cre = 0.5*log(1.+2.25*z23*xi1*ye1/massratio2) ; G4double be; if (xi <= 1.e3) be = ((2.+a6)*(1.+bet)+xi*a9)*log(1.+xii)+(a5-bet)/xi1-a9; else be = (3.-a6+a1*a7)/(2.*xi); G4double fe = (ale-cre)*be; if ( fe < 0.) fe = 0. ; G4double ymu = 4.+a6 +3.*bet*a7 ; G4double ymd = a7*(1.5+a1)*log(3.+xi)+1.-1.5*a6 ; G4double ym1 = 1.+ymu/ymd ; G4double alm_crm = log(bbb*massratio/(1.5*z23*(1.+screen*ym1))); G4double a10,bm; if ( xi >= 1.e-3) { a10 = (1.+a1)*a5 ; bm = (a7*(1.+1.5*bet)-a10*xii)*log(xi1)+xi*(a5-bet)/xi1+a10; } else { bm = (5.-a6+bet*a9)*(xi/2.); } G4double fm = alm_crm*bm; if ( fm < 0.) fm = 0. ; sum += wgi[i]*a4*(fe+fm/massratio2); } cross = -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy); return cross; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4double G4MuPairProductionModel::ComputeCrossSectionPerAtom( const G4ParticleDefinition*, G4double kineticEnergy, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy) { G4double cross = 0.0; if (kineticEnergy <= lowestKinEnergy) return cross; SetCurrentElement(Z); G4double tmax = std::min(maxEnergy, kineticEnergy); G4double cut = std::min(cutEnergy, kineticEnergy); if(cut < minPairEnergy) cut = minPairEnergy; if (cut >= tmax) return cross; cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut); if(tmax < kineticEnergy) { cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax); } return cross; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void G4MuPairProductionModel::MakeSamplingTables() { for (G4int iz=0; iz