#include <G4eBremsstrahlungRelModel.hh>
Inheritance diagram for G4eBremsstrahlungRelModel:
Definition at line 60 of file G4eBremsstrahlungRelModel.hh.
G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel | ( | const G4ParticleDefinition * | p = 0 , |
|
const G4String & | nam = "eBremLPM" | |||
) |
Definition at line 83 of file G4eBremsstrahlungRelModel.cc.
References currentZ, densityCorr, densityFactor, fParticleChange, G4Gamma::Gamma(), G4NistManager::Instance(), kinEnergy, nist, particleMass, G4VEmModel::SetAngularDistribution(), G4VEmModel::SetLowEnergyLimit(), G4VEmModel::SetLPMFlag(), theGamma, and totalEnergy.
00085 : G4VEmModel(nam), 00086 particle(0), 00087 bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.), 00088 isElectron(true), 00089 fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi), 00090 fLPMconstant(fine_structure_const*electron_mass_c2*electron_mass_c2/(4.*pi*hbarc)*0.5), 00091 fXiLPM(0), fPhiLPM(0), fGLPM(0), 00092 use_completescreening(false),isInitialised(false) 00093 { 00094 fParticleChange = 0; 00095 theGamma = G4Gamma::Gamma(); 00096 00097 lowKinEnergy = 0.1*GeV; 00098 SetLowEnergyLimit(lowKinEnergy); 00099 00100 nist = G4NistManager::Instance(); 00101 00102 SetLPMFlag(true); 00103 //SetAngularDistribution(new G4ModifiedTsai()); 00104 SetAngularDistribution(new G4DipBustGenerator()); 00105 00106 particleMass = kinEnergy = totalEnergy = currentZ = z13 = z23 = lnZ = Fel 00107 = Finel = fCoulomb = fMax = densityFactor = densityCorr = lpmEnergy 00108 = xiLPM = phiLPM = gLPM = klpm = kp = 0.0; 00109 00110 energyThresholdLPM = 1.e39; 00111 00112 InitialiseConstants(); 00113 if(p) { SetParticle(p); } 00114 }
G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel | ( | ) | [virtual] |
G4double G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom | ( | const G4ParticleDefinition * | , | |
G4double | tkin, | |||
G4double | Z, | |||
G4double | , | |||
G4double | cutEnergy, | |||
G4double | maxEnergy = DBL_MAX | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 260 of file G4eBremsstrahlungRelModel.cc.
References bremFactor, kinEnergy, particle, and SetCurrentElement().
00266 { 00267 if(!particle) { SetParticle(p); } 00268 if(kineticEnergy < lowKinEnergy) { return 0.0; } 00269 G4double cut = std::min(cutEnergy, kineticEnergy); 00270 G4double tmax = std::min(maxEnergy, kineticEnergy); 00271 00272 if(cut >= tmax) { return 0.0; } 00273 00274 SetCurrentElement(Z); 00275 00276 G4double cross = ComputeXSectionPerAtom(cut); 00277 00278 // allow partial integration 00279 if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); } 00280 00281 cross *= Z*Z*bremFactor; 00282 00283 return cross; 00284 }
G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume | ( | const G4Material * | , | |
const G4ParticleDefinition * | , | |||
G4double | kineticEnergy, | |||
G4double | cutEnergy | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 190 of file G4eBremsstrahlungRelModel.cc.
References bremFactor, currentZ, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), G4Material::GetNumberOfElements(), particle, SetCurrentElement(), G4VEmModel::SetCurrentElement(), and SetupForMaterial().
00195 { 00196 if(!particle) { SetParticle(p); } 00197 if(kineticEnergy < lowKinEnergy) { return 0.0; } 00198 G4double cut = std::min(cutEnergy, kineticEnergy); 00199 if(cut == 0.0) { return 0.0; } 00200 00201 SetupForMaterial(particle, material,kineticEnergy); 00202 00203 const G4ElementVector* theElementVector = material->GetElementVector(); 00204 const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector(); 00205 00206 G4double dedx = 0.0; 00207 00208 // loop for elements in the material 00209 for (size_t i=0; i<material->GetNumberOfElements(); i++) { 00210 00211 G4VEmModel::SetCurrentElement((*theElementVector)[i]); 00212 SetCurrentElement((*theElementVector)[i]->GetZ()); 00213 00214 dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut); 00215 } 00216 dedx *= bremFactor; 00217 00218 00219 return dedx; 00220 }
G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom | ( | G4double | gammaEnergy | ) | [protected, virtual] |
Reimplemented in G4SeltzerBergerModel.
Definition at line 428 of file G4eBremsstrahlungRelModel.cc.
References currentZ, and totalEnergy.
Referenced by SampleSecondaries().
00433 { 00434 00435 if(gammaEnergy < 0.0) { return 0.0; } 00436 00437 G4double y = gammaEnergy/totalEnergy; 00438 00439 G4double main=0.,secondTerm=0.; 00440 00441 if (use_completescreening|| currentZ<5) { 00442 // ** form factors complete screening case ** 00443 main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ ); 00444 secondTerm = (1.-y)/12.*(1.+1./currentZ); 00445 } 00446 else { 00447 // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5** 00448 G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy); 00449 G4double gg=dd/z13; 00450 G4double eps=dd/z23; 00451 G4double phi1=Phi1(gg,currentZ), phi1m2=Phi1M2(gg,currentZ); 00452 G4double psi1=Psi1(eps,currentZ), psi1m2=Psi1M2(eps,currentZ); 00453 00454 main = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ ); 00455 secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ); 00456 } 00457 G4double cross = main+secondTerm; 00458 return cross; 00459 }
void G4eBremsstrahlungRelModel::Initialise | ( | const G4ParticleDefinition * | , | |
const G4DataVector & | ||||
) | [virtual] |
Implements G4VEmModel.
Reimplemented in G4SeltzerBergerModel.
Definition at line 172 of file G4eBremsstrahlungRelModel.cc.
References currentZ, fParticleChange, G4VEmModel::GetParticleChangeForLoss(), G4VEmModel::InitialiseElementSelectors(), and G4VEmModel::LowEnergyLimit().
Referenced by G4SeltzerBergerModel::Initialise().
00174 { 00175 if(p) { SetParticle(p); } 00176 00177 lowKinEnergy = LowEnergyLimit(); 00178 00179 currentZ = 0.; 00180 00181 InitialiseElementSelectors(p, cuts); 00182 00183 if(isInitialised) { return; } 00184 fParticleChange = GetParticleChangeForLoss(); 00185 isInitialised = true; 00186 }
G4double G4eBremsstrahlungRelModel::LPMconstant | ( | ) | const [inline] |
void G4eBremsstrahlungRelModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | cutEnergy, | |||
G4double | maxEnergy | |||
) | [virtual] |
Implements G4VEmModel.
Reimplemented in G4SeltzerBergerModel.
Definition at line 463 of file G4eBremsstrahlungRelModel.cc.
References ComputeDXSectionPerAtom(), currentZ, densityCorr, densityFactor, fParticleChange, fStopAndKill, G4cout, G4endl, G4lrint(), G4UniformRand, G4VEmModel::GetAngularDistribution(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4VEmModel::GetName(), kinEnergy, particle, particleMass, G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SecondaryThreshold(), G4VEmModel::SelectRandomAtom(), SetCurrentElement(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), SetupForMaterial(), theGamma, and totalEnergy.
00469 { 00470 G4double kineticEnergy = dp->GetKineticEnergy(); 00471 if(kineticEnergy < lowKinEnergy) { return; } 00472 G4double cut = std::min(cutEnergy, kineticEnergy); 00473 G4double emax = std::min(maxEnergy, kineticEnergy); 00474 if(cut >= emax) { return; } 00475 00476 SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy); 00477 00478 const G4Element* elm = 00479 SelectRandomAtom(couple,particle,kineticEnergy,cut,emax); 00480 SetCurrentElement(elm->GetZ()); 00481 00482 kinEnergy = kineticEnergy; 00483 totalEnergy = kineticEnergy + particleMass; 00484 densityCorr = densityFactor*totalEnergy*totalEnergy; 00485 00486 //G4double fmax= fMax; 00487 G4bool highe = true; 00488 if(totalEnergy < energyThresholdLPM) { highe = false; } 00489 00490 G4double xmin = log(cut*cut + densityCorr); 00491 G4double xmax = log(emax*emax + densityCorr); 00492 G4double gammaEnergy, f, x; 00493 00494 do { 00495 x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr; 00496 if(x < 0.0) { x = 0.0; } 00497 gammaEnergy = sqrt(x); 00498 if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); } 00499 else { f = ComputeDXSectionPerAtom(gammaEnergy); } 00500 00501 if ( f > fMax ) { 00502 G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! " 00503 << f << " > " << fMax 00504 << " Egamma(MeV)= " << gammaEnergy 00505 << " Ee(MeV)= " << kineticEnergy 00506 << " " << GetName() 00507 << G4endl; 00508 } 00509 00510 } while (f < fMax*G4UniformRand()); 00511 00512 // 00513 // angles of the emitted gamma. ( Z - axis along the parent particle) 00514 // use general interface 00515 // 00516 00517 G4ThreeVector gammaDirection = 00518 GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy, 00519 G4lrint(currentZ), 00520 couple->GetMaterial()); 00521 00522 // create G4DynamicParticle object for the Gamma 00523 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection, 00524 gammaEnergy); 00525 vdp->push_back(gamma); 00526 00527 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2)); 00528 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection() 00529 - gammaEnergy*gammaDirection).unit(); 00530 00531 // energy of primary 00532 G4double finalE = kineticEnergy - gammaEnergy; 00533 00534 // stop tracking and create new secondary instead of primary 00535 if(gammaEnergy > SecondaryThreshold()) { 00536 fParticleChange->ProposeTrackStatus(fStopAndKill); 00537 fParticleChange->SetProposedKineticEnergy(0.0); 00538 G4DynamicParticle* el = 00539 new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle), 00540 direction, finalE); 00541 vdp->push_back(el); 00542 00543 // continue tracking 00544 } else { 00545 fParticleChange->SetProposedMomentumDirection(direction); 00546 fParticleChange->SetProposedKineticEnergy(finalE); 00547 } 00548 }
void G4eBremsstrahlungRelModel::SetCurrentElement | ( | const | G4double | ) | [inline, protected] |
Definition at line 185 of file G4eBremsstrahlungRelModel.hh.
References currentZ, G4VEmModel::GetCurrentElement(), G4Element::GetfCoulomb(), G4NistManager::GetLOGZ(), G4NistManager::GetZ13(), and nist.
Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), G4SeltzerBergerModel::SampleSecondaries(), and SampleSecondaries().
00186 { 00187 if(Z != currentZ) { 00188 currentZ = Z; 00189 00190 G4int iz = G4int(Z); 00191 z13 = nist->GetZ13(iz); 00192 z23 = z13*z13; 00193 lnZ = nist->GetLOGZ(iz); 00194 00195 if (iz <= 4) { 00196 Fel = Fel_light[iz]; 00197 Finel = Finel_light[iz] ; 00198 } 00199 else { 00200 Fel = facFel - lnZ/3. ; 00201 Finel = facFinel - 2.*lnZ/3. ; 00202 } 00203 00204 fCoulomb = GetCurrentElement()->GetfCoulomb(); 00205 fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.; 00206 } 00207 }
void G4eBremsstrahlungRelModel::SetLPMconstant | ( | G4double | val | ) | [inline] |
void G4eBremsstrahlungRelModel::SetupForMaterial | ( | const G4ParticleDefinition * | , | |
const G4Material * | , | |||
G4double | ||||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 145 of file G4eBremsstrahlungRelModel.cc.
References densityCorr, densityFactor, G4Material::GetElectronDensity(), G4Material::GetRadlen(), kinEnergy, G4VEmModel::LPMFlag(), particleMass, and totalEnergy.
Referenced by ComputeDEDXPerVolume(), G4SeltzerBergerModel::SampleSecondaries(), and SampleSecondaries().
00148 { 00149 densityFactor = mat->GetElectronDensity()*fMigdalConstant; 00150 lpmEnergy = mat->GetRadlen()*fLPMconstant; 00151 00152 // Threshold for LPM effect (i.e. below which LPM hidden by density effect) 00153 if (LPMFlag()) { 00154 energyThresholdLPM=sqrt(densityFactor)*lpmEnergy; 00155 } else { 00156 energyThresholdLPM=1.e39; // i.e. do not use LPM effect 00157 } 00158 // calculate threshold for density effect 00159 kinEnergy = kineticEnergy; 00160 totalEnergy = kineticEnergy + particleMass; 00161 densityCorr = densityFactor*totalEnergy*totalEnergy; 00162 00163 // define critical gamma energies (important for integration/dicing) 00164 klpm=totalEnergy*totalEnergy/lpmEnergy; 00165 kp=sqrt(densityCorr); 00166 00167 }
G4double G4eBremsstrahlungRelModel::bremFactor [protected] |
Definition at line 132 of file G4eBremsstrahlungRelModel.hh.
Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), and G4SeltzerBergerModel::ComputeDXSectionPerAtom().
G4double G4eBremsstrahlungRelModel::currentZ [protected] |
Definition at line 138 of file G4eBremsstrahlungRelModel.hh.
Referenced by ComputeDEDXPerVolume(), G4SeltzerBergerModel::ComputeDXSectionPerAtom(), ComputeDXSectionPerAtom(), G4eBremsstrahlungRelModel(), Initialise(), G4SeltzerBergerModel::SampleSecondaries(), SampleSecondaries(), and SetCurrentElement().
G4double G4eBremsstrahlungRelModel::densityCorr [protected] |
Definition at line 142 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4eBremsstrahlungRelModel(), G4SeltzerBergerModel::SampleSecondaries(), SampleSecondaries(), and SetupForMaterial().
G4double G4eBremsstrahlungRelModel::densityFactor [protected] |
Definition at line 141 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4eBremsstrahlungRelModel(), G4SeltzerBergerModel::SampleSecondaries(), SampleSecondaries(), and SetupForMaterial().
Definition at line 130 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4eBremsstrahlungRelModel(), Initialise(), G4SeltzerBergerModel::SampleSecondaries(), and SampleSecondaries().
G4bool G4eBremsstrahlungRelModel::isElectron [protected] |
Definition at line 144 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4SeltzerBergerModel::ComputeDXSectionPerAtom(), and G4SeltzerBergerModel::SampleSecondaries().
G4double G4eBremsstrahlungRelModel::kinEnergy [protected] |
Definition at line 136 of file G4eBremsstrahlungRelModel.hh.
Referenced by ComputeCrossSectionPerAtom(), G4SeltzerBergerModel::ComputeDXSectionPerAtom(), G4eBremsstrahlungRelModel(), SampleSecondaries(), and SetupForMaterial().
G4NistManager* G4eBremsstrahlungRelModel::nist [protected] |
Definition at line 127 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4eBremsstrahlungRelModel(), and SetCurrentElement().
const G4ParticleDefinition* G4eBremsstrahlungRelModel::particle [protected] |
Definition at line 128 of file G4eBremsstrahlungRelModel.hh.
Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), G4SeltzerBergerModel::SampleSecondaries(), and SampleSecondaries().
G4double G4eBremsstrahlungRelModel::particleMass [protected] |
Definition at line 135 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4SeltzerBergerModel::ComputeDXSectionPerAtom(), G4eBremsstrahlungRelModel(), G4SeltzerBergerModel::SampleSecondaries(), SampleSecondaries(), and SetupForMaterial().
G4ParticleDefinition* G4eBremsstrahlungRelModel::theGamma [protected] |
Definition at line 129 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4eBremsstrahlungRelModel(), G4SeltzerBergerModel::SampleSecondaries(), and SampleSecondaries().
G4double G4eBremsstrahlungRelModel::totalEnergy [protected] |
Definition at line 137 of file G4eBremsstrahlungRelModel.hh.
Referenced by G4SeltzerBergerModel::ComputeDXSectionPerAtom(), ComputeDXSectionPerAtom(), G4eBremsstrahlungRelModel(), G4SeltzerBergerModel::SampleSecondaries(), SampleSecondaries(), and SetupForMaterial().