Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeAnnihilationModel Class Reference

#include <G4PenelopeAnnihilationModel.hh>

+ Inheritance diagram for G4PenelopeAnnihilationModel:

Public Member Functions

 G4PenelopeAnnihilationModel (const G4ParticleDefinition *p=0, const G4String &processName="PenAnnih")
 
virtual ~G4PenelopeAnnihilationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
const G4ParticleDefinitionfParticle
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 55 of file G4PenelopeAnnihilationModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeAnnihilationModel()

G4PenelopeAnnihilationModel::G4PenelopeAnnihilationModel ( const G4ParticleDefinition p = 0,
const G4String processName = "PenAnnih" 
)

Definition at line 52 of file G4PenelopeAnnihilationModel.cc.

54 :G4VEmModel(nam),fParticleChange(0),fParticle(0),isInitialised(false)
55{
56 fIntrinsicLowEnergyLimit = 0.0;
57 fIntrinsicHighEnergyLimit = 100.0*GeV;
58 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
59 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
60
61 if (part)
62 SetParticle(part);
63
64 //Calculate variable that will be used later on
65 fPielr2 = pi*classic_electr_radius*classic_electr_radius;
66
67 verboseLevel= 0;
68 // Verbosity scale:
69 // 0 = nothing
70 // 1 = warning for energy non-conservation
71 // 2 = details of energy budget
72 // 3 = calculation of cross sections, file openings, sampling of atoms
73 // 4 = entering in methods
74
75}
const G4ParticleDefinition * fParticle
G4ParticleChangeForGamma * fParticleChange
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
const G4double pi

◆ ~G4PenelopeAnnihilationModel()

G4PenelopeAnnihilationModel::~G4PenelopeAnnihilationModel ( )
virtual

Definition at line 79 of file G4PenelopeAnnihilationModel.cc.

80{;}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PenelopeAnnihilationModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 132 of file G4PenelopeAnnihilationModel.cc.

137{
138 if (verboseLevel > 3)
139 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopeAnnihilationModel" <<
140 G4endl;
141
142 G4double cs = Z*ComputeCrossSectionPerElectron(energy);
143
144 if (verboseLevel > 2)
145 G4cout << "Annihilation cross Section at " << energy/keV << " keV for Z=" << Z <<
146 " = " << cs/barn << " barn" << G4endl;
147 return cs;
148}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double energy(const ThreeVector &p, const G4double m)

◆ GetVerbosityLevel()

G4int G4PenelopeAnnihilationModel::GetVerbosityLevel ( )
inline

Definition at line 82 of file G4PenelopeAnnihilationModel.hh.

82{return verboseLevel;};

◆ Initialise()

void G4PenelopeAnnihilationModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 84 of file G4PenelopeAnnihilationModel.cc.

86{
87 if (verboseLevel > 3)
88 G4cout << "Calling G4PenelopeAnnihilationModel::Initialise()" << G4endl;
89 SetParticle(part);
90
91 if (IsMaster() && part == fParticle)
92 {
93
94 if(verboseLevel > 0) {
95 G4cout << "Penelope Annihilation model is initialized " << G4endl
96 << "Energy range: "
97 << LowEnergyLimit() / keV << " keV - "
98 << HighEnergyLimit() / GeV << " GeV"
99 << G4endl;
100 }
101 }
102
103 if(isInitialised) return;
105 isInitialised = true;
106}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ InitialiseLocal()

void G4PenelopeAnnihilationModel::InitialiseLocal ( const G4ParticleDefinition part,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 109 of file G4PenelopeAnnihilationModel.cc.

111{
112 if (verboseLevel > 3)
113 G4cout << "Calling G4PenelopeAnnihilationModel::InitialiseLocal()" << G4endl;
114
115 //
116 //Check that particle matches: one might have multiple master models (e.g.
117 //for e+ and e-).
118 //
119 if (part == fParticle)
120 {
121 //Get the const table pointers from the master to the workers
122 const G4PenelopeAnnihilationModel* theModel =
123 static_cast<G4PenelopeAnnihilationModel*> (masterModel);
124
125 //Same verbosity for all workers, as the master
126 verboseLevel = theModel->verboseLevel;
127 }
128}

◆ SampleSecondaries()

void G4PenelopeAnnihilationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicPositron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 152 of file G4PenelopeAnnihilationModel.cc.

157{
158 //
159 // Penelope model to sample final state for positron annihilation.
160 // Target eletrons are assumed to be free and at rest. Binding effects enabling
161 // one-photon annihilation are neglected.
162 // For annihilation at rest, two back-to-back photons are emitted, having energy of 511 keV
163 // and isotropic angular distribution.
164 // For annihilation in flight, it is used the theory from
165 // W. Heitler, The quantum theory of radiation, Oxford University Press (1954)
166 // The two photons can have different energy. The efficiency of the sampling algorithm
167 // of the photon energy from the dSigma/dE distribution is practically 100% for
168 // positrons of kinetic energy < 10 keV. It reaches a minimum (about 80%) at energy
169 // of about 10 MeV.
170 // The angle theta is kinematically linked to the photon energy, to ensure momentum
171 // conservation. The angle phi is sampled isotropically for the first gamma.
172 //
173 if (verboseLevel > 3)
174 G4cout << "Calling SamplingSecondaries() of G4PenelopeAnnihilationModel" << G4endl;
175
176 G4double kineticEnergy = aDynamicPositron->GetKineticEnergy();
177
178 // kill primary
181
182 if (kineticEnergy == 0.0)
183 {
184 //Old AtRestDoIt
185 G4double cosTheta = -1.0+2.0*G4UniformRand();
186 G4double sinTheta = std::sqrt(1.0-cosTheta*cosTheta);
187 G4double phi = twopi*G4UniformRand();
188 G4ThreeVector direction (sinTheta*std::cos(phi),sinTheta*std::sin(phi),cosTheta);
190 direction, electron_mass_c2);
192 -direction, electron_mass_c2);
193
194 fvect->push_back(firstGamma);
195 fvect->push_back(secondGamma);
196 return;
197 }
198
199 //This is the "PostStep" case (annihilation in flight)
200 G4ParticleMomentum positronDirection =
201 aDynamicPositron->GetMomentumDirection();
202 G4double gamma = 1.0 + std::max(kineticEnergy,1.0*eV)/electron_mass_c2;
203 G4double gamma21 = std::sqrt(gamma*gamma-1);
204 G4double ani = 1.0+gamma;
205 G4double chimin = 1.0/(ani+gamma21);
206 G4double rchi = (1.0-chimin)/chimin;
207 G4double gt0 = ani*ani-2.0;
208 G4double test=0.0;
209 G4double epsilon = 0;
210 do{
211 epsilon = chimin*std::pow(rchi,G4UniformRand());
212 G4double reject = ani*ani*(1.0-epsilon)+2.0*gamma-(1.0/epsilon);
213 test = G4UniformRand()*gt0-reject;
214 }while(test>0);
215
216 G4double totalAvailableEnergy = kineticEnergy + 2.0*electron_mass_c2;
217 G4double photon1Energy = epsilon*totalAvailableEnergy;
218 G4double photon2Energy = (1.0-epsilon)*totalAvailableEnergy;
219 G4double cosTheta1 = (ani-1.0/epsilon)/gamma21;
220 G4double cosTheta2 = (ani-1.0/(1.0-epsilon))/gamma21;
221
222 //G4double localEnergyDeposit = 0.;
223
224 G4double sinTheta1 = std::sqrt(1.-cosTheta1*cosTheta1);
225 G4double phi1 = twopi * G4UniformRand();
226 G4double dirx1 = sinTheta1 * std::cos(phi1);
227 G4double diry1 = sinTheta1 * std::sin(phi1);
228 G4double dirz1 = cosTheta1;
229
230 G4double sinTheta2 = std::sqrt(1.-cosTheta2*cosTheta2);
231 G4double phi2 = phi1+pi;
232 G4double dirx2 = sinTheta2 * std::cos(phi2);
233 G4double diry2 = sinTheta2 * std::sin(phi2);
234 G4double dirz2 = cosTheta2;
235
236 G4ThreeVector photon1Direction (dirx1,diry1,dirz1);
237 photon1Direction.rotateUz(positronDirection);
238 // create G4DynamicParticle object for the particle1
240 photon1Direction,
241 photon1Energy);
242 fvect->push_back(aParticle1);
243
244 G4ThreeVector photon2Direction(dirx2,diry2,dirz2);
245 photon2Direction.rotateUz(positronDirection);
246 // create G4DynamicParticle object for the particle2
248 photon2Direction,
249 photon2Energy);
250 fvect->push_back(aParticle2);
251
252 if (verboseLevel > 1)
253 {
254 G4cout << "-----------------------------------------------------------" << G4endl;
255 G4cout << "Energy balance from G4PenelopeAnnihilation" << G4endl;
256 G4cout << "Kinetic positron energy: " << kineticEnergy/keV << " keV" << G4endl;
257 G4cout << "Total available energy: " << totalAvailableEnergy/keV << " keV " << G4endl;
258 G4cout << "-----------------------------------------------------------" << G4endl;
259 G4cout << "Photon energy 1: " << photon1Energy/keV << " keV" << G4endl;
260 G4cout << "Photon energy 2: " << photon2Energy/keV << " keV" << G4endl;
261 G4cout << "Total final state: " << (photon1Energy+photon2Energy)/keV <<
262 " keV" << G4endl;
263 G4cout << "-----------------------------------------------------------" << G4endl;
264 }
265 if (verboseLevel > 0)
266 {
267 G4double energyDiff = std::fabs(totalAvailableEnergy-photon1Energy-photon2Energy);
268 if (energyDiff > 0.05*keV)
269 G4cout << "Warning from G4PenelopeAnnihilation: problem with energy conservation: " <<
270 (photon1Energy+photon2Energy)/keV <<
271 " keV (final) vs. " <<
272 totalAvailableEnergy/keV << " keV (initial)" << G4endl;
273 }
274 return;
275}
double epsilon(double density, double temperature)
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)

◆ SetVerbosityLevel()

void G4PenelopeAnnihilationModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 81 of file G4PenelopeAnnihilationModel.hh.

81{verboseLevel = lev;};

Member Data Documentation

◆ fParticle

const G4ParticleDefinition* G4PenelopeAnnihilationModel::fParticle
protected

Definition at line 86 of file G4PenelopeAnnihilationModel.hh.

Referenced by Initialise(), and InitialiseLocal().

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopeAnnihilationModel::fParticleChange
protected

Definition at line 85 of file G4PenelopeAnnihilationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


The documentation for this class was generated from the following files: