#include <G4IonCoulombCrossSection.hh>
Public Member Functions | |
G4IonCoulombCrossSection () | |
virtual | ~G4IonCoulombCrossSection () |
void | Initialise (const G4ParticleDefinition *, G4double cosThetaLim) |
G4double | NuclearCrossSection () |
G4double | SampleCosineTheta () |
void | SetupParticle (const G4ParticleDefinition *) |
void | SetupKinematic (G4double kinEnergy, G4double cut, G4int iz) |
void | SetupTarget (G4double Z, G4double kinEnergy, G4int heavycorr) |
G4double | GetMomentum2 () |
Protected Attributes | |
G4double | coeff |
G4double | cosThetaMin |
G4double | cosThetaMax |
G4double | cosTetMinNuc |
G4double | cosTetMaxNuc |
G4double | nucXSection |
G4double | ecut |
G4double | etag |
const G4ParticleDefinition * | particle |
G4double | chargeSquare |
G4double | spin |
G4double | mass |
G4double | tkinLab |
G4double | momLab2 |
G4double | invbetaLab2 |
G4double | tkin |
G4double | mom2 |
G4double | invbeta2 |
G4double | targetZ |
G4double | targetMass |
G4double | screenZ |
Definition at line 76 of file G4IonCoulombCrossSection.hh.
G4IonCoulombCrossSection::G4IonCoulombCrossSection | ( | ) |
Definition at line 67 of file G4IonCoulombCrossSection.cc.
References chargeSquare, coeff, cosTetMaxNuc, cosTetMinNuc, ecut, etag, G4NistManager::Instance(), invbeta2, invbetaLab2, mass, mom2, momLab2, nucXSection, particle, G4Proton::Proton(), screenZ, spin, targetMass, targetZ, tkin, and tkinLab.
00067 : 00068 cosThetaMin(1.0), 00069 cosThetaMax(-1.0), 00070 alpha2(fine_structure_const*fine_structure_const) 00071 { 00072 fNistManager = G4NistManager::Instance(); 00073 theProton = G4Proton::Proton(); 00074 particle=0; 00075 00076 G4double p0 = electron_mass_c2*classic_electr_radius; 00077 coeff = twopi*p0*p0; 00078 00079 cosTetMinNuc=0; 00080 cosTetMaxNuc=0; 00081 nucXSection =0; 00082 00083 00084 chargeSquare = spin = mass =0; 00085 tkinLab = momLab2 = invbetaLab2=0; 00086 tkin = mom2 = invbeta2=0; 00087 00088 targetZ = targetMass = screenZ =0; 00089 ScreenRSquare=0.; 00090 00091 etag = ecut = 0.0; 00092 00093 } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4IonCoulombCrossSection::~G4IonCoulombCrossSection | ( | ) | [virtual] |
G4double G4IonCoulombCrossSection::GetMomentum2 | ( | ) | [inline] |
Definition at line 177 of file G4IonCoulombCrossSection.hh.
References mom2.
Referenced by G4IonCoulombScatteringModel::SampleSecondaries().
00177 { 00178 return mom2; 00179 }
void G4IonCoulombCrossSection::Initialise | ( | const G4ParticleDefinition * | , | |
G4double | cosThetaLim | |||
) |
Definition at line 100 of file G4IonCoulombCrossSection.cc.
References cosThetaMin, DBL_MAX, DBL_MIN, ecut, etag, mom2, nucXSection, particle, SetupParticle(), targetZ, and tkin.
Referenced by G4IonCoulombScatteringModel::Initialise().
00102 { 00103 00104 SetupParticle(p); 00105 nucXSection = 0.0; 00106 tkin = targetZ = mom2 = DBL_MIN; 00107 ecut = etag = DBL_MAX; 00108 particle = p; 00109 00110 00111 cosThetaMin = CosThetaLim; 00112 00113 }
G4double G4IonCoulombCrossSection::NuclearCrossSection | ( | ) |
Definition at line 211 of file G4IonCoulombCrossSection.cc.
References chargeSquare, coeff, cosTetMaxNuc, cosTetMinNuc, invbeta2, mom2, nucXSection, screenZ, and targetZ.
Referenced by G4IonCoulombScatteringModel::ComputeCrossSectionPerAtom().
00212 { 00213 // This method needs initialisation before be called 00214 00215 // scattering with target nucleus 00216 G4double fac = coeff*targetZ*targetZ*chargeSquare*invbeta2/mom2; 00217 00218 nucXSection = 0.0; 00219 00220 G4double x = 1.0 - cosTetMinNuc; 00221 G4double x1 = x + screenZ; 00222 00223 // scattering with nucleus 00224 if(cosTetMaxNuc < cosTetMinNuc) { 00225 00226 nucXSection =fac*(cosTetMinNuc - cosTetMaxNuc)/ 00227 (x1*(1.0 - cosTetMaxNuc + screenZ)); 00228 00229 } 00230 00231 return nucXSection; 00232 }
G4double G4IonCoulombCrossSection::SampleCosineTheta | ( | ) |
Definition at line 236 of file G4IonCoulombCrossSection.cc.
References cosTetMaxNuc, cosTetMinNuc, G4UniformRand, and screenZ.
Referenced by G4IonCoulombScatteringModel::SampleSecondaries().
00237 { 00238 00239 if(cosTetMaxNuc >= cosTetMinNuc) return 0.0; 00240 00241 G4double x1 = 1. - cosTetMinNuc + screenZ; 00242 G4double x2 = 1. - cosTetMaxNuc + screenZ; 00243 G4double dx = cosTetMinNuc - cosTetMaxNuc; 00244 G4double /*grej,*/ z1; 00245 00246 z1 = x1*x2/(x1 + G4UniformRand()*dx) - screenZ; 00247 //grej = 1.0/(1.0 + z1); 00248 return z1; 00249 }
Definition at line 116 of file G4IonCoulombCrossSection.cc.
References cosTetMaxNuc, cosTetMinNuc, cosThetaMax, cosThetaMin, ecut, G4NistManager::GetAtomicMassAmu(), invbeta2, invbetaLab2, mass, mom2, momLab2, targetMass, tkin, and tkinLab.
Referenced by G4IonCoulombScatteringModel::ComputeCrossSectionPerAtom().
00118 { 00119 if(ekin != tkinLab || ecut != cut) { 00120 00121 // lab 00122 tkinLab = ekin; 00123 momLab2 = tkinLab*(tkinLab + 2.0*mass); 00124 invbetaLab2 = 1.0 + mass*mass/momLab2; 00125 00126 G4double etot = tkinLab + mass; 00127 G4double ptot = sqrt(momLab2); 00128 G4double m12 = mass*mass; 00129 00130 targetMass=fNistManager->GetAtomicMassAmu(iz)*amu_c2; 00131 00132 // relativistic reduced mass from publucation 00133 // A.P. Martynenko, R.N. Faustov, Teoret. mat. Fiz. 64 (1985) 179 00134 00135 //incident particle & target nucleus 00136 G4double Ecm=sqrt(m12 + targetMass*targetMass + 2.0*etot*targetMass); 00137 G4double mu_rel=mass*targetMass/Ecm; 00138 G4double momCM= ptot*targetMass/Ecm; 00139 // relative system 00140 mom2 = momCM*momCM; 00141 invbeta2 = 1.0 + mu_rel*mu_rel/mom2; 00142 tkin = momCM*sqrt(invbeta2) - mu_rel;//Ekin of mu_rel 00143 //......................................................... 00144 00145 cosTetMinNuc = cosThetaMin; 00146 cosTetMaxNuc = cosThetaMax; 00147 00148 } 00149 00150 }
void G4IonCoulombCrossSection::SetupParticle | ( | const G4ParticleDefinition * | ) | [inline] |
Definition at line 164 of file G4IonCoulombCrossSection.hh.
References chargeSquare, G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4ParticleDefinition::GetPDGSpin(), mass, particle, spin, and tkin.
Referenced by Initialise(), and G4IonCoulombScatteringModel::SetupParticle().
00165 { 00166 particle = p; 00167 mass = particle->GetPDGMass(); 00168 spin = particle->GetPDGSpin(); 00169 if(0.0 != spin) { spin = 0.5; } 00170 G4double q = std::fabs(particle->GetPDGCharge()/CLHEP::eplus); 00171 chargeSquare = q*q; 00172 tkin = 0.0; 00173 }
Definition at line 154 of file G4IonCoulombCrossSection.cc.
References chargeSquare, cosTetMaxNuc, etag, invbeta2, mom2, particle, screenZ, and targetZ.
Referenced by G4IonCoulombScatteringModel::ComputeCrossSectionPerAtom().
00155 { 00156 if(Z != targetZ || e != etag) { 00157 etag = e; 00158 targetZ = Z; 00159 G4int iz= G4int(Z); 00160 00161 SetScreenRSquare(iz); 00162 screenZ =0; 00163 screenZ = ScreenRSquare/mom2; 00164 00165 // G4cout<< "heavycorr "<<heavycorr<<G4endl; 00166 00167 if(heavycorr!=0 && particle != theProton){ 00168 G4double corr=5.*twopi*Z*std::sqrt(chargeSquare*alpha2); 00169 corr=std::pow(corr,0.12); 00170 screenZ *=(1.13 + corr*3.76*Z*Z*chargeSquare*invbeta2*alpha2)/2.; 00171 // G4cout<<" heavycorr Z e corr....2As "<< heavycorr << "\t" 00172 // <<Z <<"\t"<<e/MeV <<"\t"<<screenZ<<G4endl; 00173 00174 }else{ screenZ *=(1.13 + 3.76*Z*Z*chargeSquare*invbeta2*alpha2)/2.; 00175 // G4cout<<" heavycorr Z e....2As "<< heavycorr << "\t" 00176 // <<Z <<"\t"<< e/MeV <<"\t" <<screenZ<<G4endl; 00177 } 00178 00179 if(1 == iz && particle == theProton && cosTetMaxNuc < 0.0) { 00180 cosTetMaxNuc = 0.0; 00181 } 00182 00183 } 00184 }
G4double G4IonCoulombCrossSection::chargeSquare [protected] |
Definition at line 134 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), NuclearCrossSection(), SetupParticle(), and SetupTarget().
G4double G4IonCoulombCrossSection::coeff [protected] |
Definition at line 113 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), and NuclearCrossSection().
G4double G4IonCoulombCrossSection::cosTetMaxNuc [protected] |
Definition at line 121 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), NuclearCrossSection(), SampleCosineTheta(), SetupKinematic(), and SetupTarget().
G4double G4IonCoulombCrossSection::cosTetMinNuc [protected] |
Definition at line 120 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), NuclearCrossSection(), SampleCosineTheta(), and SetupKinematic().
G4double G4IonCoulombCrossSection::cosThetaMax [protected] |
G4double G4IonCoulombCrossSection::cosThetaMin [protected] |
Definition at line 116 of file G4IonCoulombCrossSection.hh.
Referenced by Initialise(), and SetupKinematic().
G4double G4IonCoulombCrossSection::ecut [protected] |
Definition at line 128 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), Initialise(), and SetupKinematic().
G4double G4IonCoulombCrossSection::etag [protected] |
Definition at line 129 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), Initialise(), and SetupTarget().
G4double G4IonCoulombCrossSection::invbeta2 [protected] |
Definition at line 146 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), NuclearCrossSection(), SetupKinematic(), and SetupTarget().
G4double G4IonCoulombCrossSection::invbetaLab2 [protected] |
Definition at line 141 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), and SetupKinematic().
G4double G4IonCoulombCrossSection::mass [protected] |
Definition at line 136 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), SetupKinematic(), and SetupParticle().
G4double G4IonCoulombCrossSection::mom2 [protected] |
Definition at line 145 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), GetMomentum2(), Initialise(), NuclearCrossSection(), SetupKinematic(), and SetupTarget().
G4double G4IonCoulombCrossSection::momLab2 [protected] |
Definition at line 140 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), and SetupKinematic().
G4double G4IonCoulombCrossSection::nucXSection [protected] |
Definition at line 125 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), Initialise(), and NuclearCrossSection().
const G4ParticleDefinition* G4IonCoulombCrossSection::particle [protected] |
Definition at line 132 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), Initialise(), SetupParticle(), and SetupTarget().
G4double G4IonCoulombCrossSection::screenZ [protected] |
Definition at line 151 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), NuclearCrossSection(), SampleCosineTheta(), and SetupTarget().
G4double G4IonCoulombCrossSection::spin [protected] |
Definition at line 135 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), and SetupParticle().
G4double G4IonCoulombCrossSection::targetMass [protected] |
Definition at line 150 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), and SetupKinematic().
G4double G4IonCoulombCrossSection::targetZ [protected] |
Definition at line 149 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), Initialise(), NuclearCrossSection(), and SetupTarget().
G4double G4IonCoulombCrossSection::tkin [protected] |
Definition at line 144 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), Initialise(), SetupKinematic(), and SetupParticle().
G4double G4IonCoulombCrossSection::tkinLab [protected] |
Definition at line 139 of file G4IonCoulombCrossSection.hh.
Referenced by G4IonCoulombCrossSection(), and SetupKinematic().