Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PenelopeAnnihilationModel Class Reference

#include <G4PenelopeAnnihilationModel.hh>

Inheritance diagram for G4PenelopeAnnihilationModel:
Collaboration 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 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 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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, 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 *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Attributes

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

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 56 of file G4PenelopeAnnihilationModel.hh.

Constructor & Destructor Documentation

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

Definition at line 53 of file G4PenelopeAnnihilationModel.cc.

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

Here is the call graph for this function:

G4PenelopeAnnihilationModel::~G4PenelopeAnnihilationModel ( )
virtual

Definition at line 80 of file G4PenelopeAnnihilationModel.cc.

81 {;}

Member Function Documentation

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 133 of file G4PenelopeAnnihilationModel.cc.

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

Definition at line 83 of file G4PenelopeAnnihilationModel.hh.

83 {return verboseLevel;};
void G4PenelopeAnnihilationModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 85 of file G4PenelopeAnnihilationModel.cc.

87 {
88  if (verboseLevel > 3)
89  G4cout << "Calling G4PenelopeAnnihilationModel::Initialise()" << G4endl;
90  SetParticle(part);
91 
92  if (IsMaster() && part == fParticle)
93  {
94 
95  if(verboseLevel > 0) {
96  G4cout << "Penelope Annihilation model is initialized " << G4endl
97  << "Energy range: "
98  << LowEnergyLimit() / keV << " keV - "
99  << HighEnergyLimit() / GeV << " GeV"
100  << G4endl;
101  }
102  }
103 
104  if(isInitialised) return;
106  isInitialised = true;
107 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4ParticleDefinition * fParticle
G4GLOB_DLL std::ostream G4cout
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * fParticleChange
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 110 of file G4PenelopeAnnihilationModel.cc.

112 {
113  if (verboseLevel > 3)
114  G4cout << "Calling G4PenelopeAnnihilationModel::InitialiseLocal()" << G4endl;
115 
116  //
117  //Check that particle matches: one might have multiple master models (e.g.
118  //for e+ and e-).
119  //
120  if (part == fParticle)
121  {
122  //Get the const table pointers from the master to the workers
123  const G4PenelopeAnnihilationModel* theModel =
124  static_cast<G4PenelopeAnnihilationModel*> (masterModel);
125 
126  //Same verbosity for all workers, as the master
127  verboseLevel = theModel->verboseLevel;
128  }
129 }
const G4ParticleDefinition * fParticle
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4PenelopeAnnihilationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicPositron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 153 of file G4PenelopeAnnihilationModel.cc.

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

Here is the call graph for this function:

void G4PenelopeAnnihilationModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 82 of file G4PenelopeAnnihilationModel.hh.

82 {verboseLevel = lev;};

Member Data Documentation

const G4ParticleDefinition* G4PenelopeAnnihilationModel::fParticle
protected

Definition at line 87 of file G4PenelopeAnnihilationModel.hh.

G4ParticleChangeForGamma* G4PenelopeAnnihilationModel::fParticleChange
protected

Definition at line 83 of file G4PenelopeAnnihilationModel.hh.


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