#include <G4eplusAnnihilation.hh>
Inheritance diagram for G4eplusAnnihilation:
Public Member Functions | |
G4eplusAnnihilation (const G4String &name="annihil") | |
virtual | ~G4eplusAnnihilation () |
virtual G4bool | IsApplicable (const G4ParticleDefinition &p) |
virtual G4VParticleChange * | AtRestDoIt (const G4Track &track, const G4Step &stepData) |
virtual G4double | AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) |
virtual void | PrintInfo () |
Protected Member Functions | |
virtual void | InitialiseProcess (const G4ParticleDefinition *) |
Definition at line 67 of file G4eplusAnnihilation.hh.
G4eplusAnnihilation::G4eplusAnnihilation | ( | const G4String & | name = "annihil" |
) |
Definition at line 68 of file G4eplusAnnihilation.cc.
References G4VProcess::enableAtRestDoIt, fAnnihilation, G4Gamma::Gamma(), G4VEmProcess::SetBuildTableFlag(), G4VEmProcess::SetIntegral(), G4VProcess::SetProcessSubType(), G4VEmProcess::SetSecondaryParticle(), and G4VEmProcess::SetStartFromNullFlag().
00069 : G4VEmProcess(name), isInitialised(false) 00070 { 00071 theGamma = G4Gamma::Gamma(); 00072 SetIntegral(true); 00073 SetBuildTableFlag(false); 00074 SetStartFromNullFlag(false); 00075 SetSecondaryParticle(theGamma); 00076 SetProcessSubType(fAnnihilation); 00077 enableAtRestDoIt = true; 00078 }
G4eplusAnnihilation::~G4eplusAnnihilation | ( | ) | [virtual] |
G4VParticleChange * G4eplusAnnihilation::AtRestDoIt | ( | const G4Track & | track, | |
const G4Step & | stepData | |||
) | [virtual] |
Reimplemented from G4VDiscreteProcess.
Definition at line 121 of file G4eplusAnnihilation.cc.
References G4ParticleChangeForGamma::AddSecondary(), G4VEmProcess::fParticleChange, fStopAndKill, G4UniformRand, G4Track::GetWeight(), G4ParticleChangeForGamma::InitializeForPostStep(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::ProposeWeight(), G4VParticleChange::SetNumberOfSecondaries(), and G4DynamicParticle::SetPolarization().
00129 : Effects due to binding of atomic electrons are negliged. 00130 { 00131 fParticleChange.InitializeForPostStep(aTrack); 00132 00133 G4double cosTeta = 2.*G4UniformRand()-1.; 00134 G4double sinTeta = sqrt((1.-cosTeta)*(1.0 + cosTeta)); 00135 G4double phi = twopi * G4UniformRand(); 00136 G4ThreeVector dir(sinTeta*cos(phi), sinTeta*sin(phi), cosTeta); 00137 phi = twopi * G4UniformRand(); 00138 G4ThreeVector pol(cos(phi), sin(phi), 0.0); 00139 pol.rotateUz(dir); 00140 00141 fParticleChange.ProposeWeight(aTrack.GetWeight()); 00142 // add gammas 00143 fParticleChange.SetNumberOfSecondaries(2); 00144 G4DynamicParticle* dp = 00145 new G4DynamicParticle(theGamma, dir, electron_mass_c2); 00146 dp->SetPolarization(pol.x(),pol.y(),pol.z()); 00147 fParticleChange.AddSecondary(dp); 00148 00149 dp = new G4DynamicParticle(theGamma,-dir, electron_mass_c2); 00150 dp->SetPolarization(-pol.x(),-pol.y(),-pol.z()); 00151 fParticleChange.AddSecondary(dp); 00152 00153 // Kill the incident positron 00154 // 00155 fParticleChange.ProposeTrackStatus(fStopAndKill); 00156 return &fParticleChange; 00157 }
G4double G4eplusAnnihilation::AtRestGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4ForceCondition * | condition | |||
) | [virtual] |
Reimplemented from G4VDiscreteProcess.
Definition at line 94 of file G4eplusAnnihilation.cc.
References NotForced.
void G4eplusAnnihilation::InitialiseProcess | ( | const G4ParticleDefinition * | ) | [protected, virtual] |
Implements G4VEmProcess.
Definition at line 103 of file G4eplusAnnihilation.cc.
References G4VEmProcess::AddEmModel(), G4VEmProcess::EmModel(), G4VEmProcess::MaxKinEnergy(), G4VEmProcess::MinKinEnergy(), G4VEmProcess::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().
00104 { 00105 if(!isInitialised) { 00106 isInitialised = true; 00107 if(!EmModel(1)) { SetEmModel(new G4eeToTwoGammaModel(),1); } 00108 EmModel(1)->SetLowEnergyLimit(MinKinEnergy()); 00109 EmModel(1)->SetHighEnergyLimit(MaxKinEnergy()); 00110 AddEmModel(1, EmModel(1)); 00111 } 00112 }
G4bool G4eplusAnnihilation::IsApplicable | ( | const G4ParticleDefinition & | p | ) | [virtual] |
Implements G4VEmProcess.
Definition at line 87 of file G4eplusAnnihilation.cc.
References G4Positron::Positron().
00088 { 00089 return (&p == G4Positron::Positron()); 00090 }
void G4eplusAnnihilation::PrintInfo | ( | ) | [virtual] |