Geant4 11.1.1
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=nullptr, const G4String &processName="PenAnnih")
 
virtual ~G4PenelopeAnnihilationModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
G4PenelopeAnnihilationModeloperator= (const G4PenelopeAnnihilationModel &right)=delete
 
 G4PenelopeAnnihilationModel (const G4PenelopeAnnihilationModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 *)
 
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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

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() [1/2]

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

Definition at line 52 of file G4PenelopeAnnihilationModel.cc.

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

◆ ~G4PenelopeAnnihilationModel()

G4PenelopeAnnihilationModel::~G4PenelopeAnnihilationModel ( )
virtual

Definition at line 77 of file G4PenelopeAnnihilationModel.cc.

78{;}

◆ G4PenelopeAnnihilationModel() [2/2]

G4PenelopeAnnihilationModel::G4PenelopeAnnihilationModel ( const G4PenelopeAnnihilationModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 130 of file G4PenelopeAnnihilationModel.cc.

135{
136 if (fVerboseLevel > 3)
137 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopeAnnihilationModel" <<
138 G4endl;
139
140 G4double cs = Z*ComputeCrossSectionPerElectron(energy);
141
142 if (fVerboseLevel > 2)
143 G4cout << "Annihilation cross Section at " << energy/keV << " keV for Z=" << Z <<
144 " = " << cs/barn << " barn" << G4endl;
145 return cs;
146}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
#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 80 of file G4PenelopeAnnihilationModel.hh.

80{return fVerboseLevel;};

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 82 of file G4PenelopeAnnihilationModel.cc.

84{
85 if (fVerboseLevel > 3)
86 G4cout << "Calling G4PenelopeAnnihilationModel::Initialise()" << G4endl;
87 SetParticle(part);
88
89 if (IsMaster() && part == fParticle)
90 {
91
92 if(fVerboseLevel > 0) {
93 G4cout << "Penelope Annihilation model is initialized " << G4endl
94 << "Energy range: "
95 << LowEnergyLimit() / keV << " keV - "
96 << HighEnergyLimit() / GeV << " GeV"
97 << G4endl;
98 }
99 }
100
101 if(fIsInitialised) return;
103 fIsInitialised = true;
104}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 107 of file G4PenelopeAnnihilationModel.cc.

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

◆ operator=()

G4PenelopeAnnihilationModel & G4PenelopeAnnihilationModel::operator= ( const G4PenelopeAnnihilationModel right)
delete

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 150 of file G4PenelopeAnnihilationModel.cc.

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

79{fVerboseLevel = lev;};

Member Data Documentation

◆ fParticle

const G4ParticleDefinition* G4PenelopeAnnihilationModel::fParticle
protected

Definition at line 87 of file G4PenelopeAnnihilationModel.hh.

Referenced by Initialise(), and InitialiseLocal().

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopeAnnihilationModel::fParticleChange
protected

Definition at line 86 of file G4PenelopeAnnihilationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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