Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eplusAnnihilation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4eplusAnnihilation
34//
35// Author: Vladimir Ivanchenko on base of Michel Maire code
36//
37// Creation date: 02.08.2004
38//
39// Modifications:
40// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
41// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
42// 03-05-05 suppress Integral option (mma)
43// 04-05-05 Make class to be default (V.Ivanchenko)
44// 25-01-06 remove cut dependance in AtRestDoIt (mma)
45// 09-08-06 add SetModel(G4VEmModel*) (mma)
46// 12-09-06, move SetModel(G4VEmModel*) in G4VEmProcess (mma)
47// 30-05-12 propagate parent weight to secondaries (D. Sawkey)
48//
49
50//
51// -------------------------------------------------------------------
52//
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
59#include "G4Gamma.hh"
60#include "G4PhysicsVector.hh"
61#include "G4PhysicsLogVector.hh"
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66using namespace std;
67
69 : G4VEmProcess(name), isInitialised(false)
70{
71 theGamma = G4Gamma::Gamma();
72 SetIntegral(true);
73 SetBuildTableFlag(false);
75 SetSecondaryParticle(theGamma);
77 enableAtRestDoIt = true;
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
83{}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
88{
89 return (&p == G4Positron::Positron());
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
96{
98 return 0.0;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102
104{
105 if(!isInitialised) {
106 isInitialised = true;
107 if(!EmModel(1)) { SetEmModel(new G4eeToTwoGammaModel(),1); }
110 AddEmModel(1, EmModel(1));
111 }
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115
117{}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
122 const G4Step& )
123//
124// Performs the e+ e- annihilation when both particles are assumed at rest.
125// It generates two back to back photons with energy = electron_mass.
126// The angular distribution is isotropic.
127// GEANT4 internal units
128//
129// Note : Effects due to binding of atomic electrons are negliged.
130{
132
133 G4double cosTeta = 2.*G4UniformRand()-1.;
134 G4double sinTeta = sqrt((1.-cosTeta)*(1.0 + cosTeta));
135 G4double phi = twopi * G4UniformRand();
136 G4ThreeVector dir(sinTeta*cos(phi), sinTeta*sin(phi), cosTeta);
137 phi = twopi * G4UniformRand();
138 G4ThreeVector pol(cos(phi), sin(phi), 0.0);
139 pol.rotateUz(dir);
140
142 // add gammas
144 G4DynamicParticle* dp =
145 new G4DynamicParticle(theGamma, dir, electron_mass_c2);
146 dp->SetPolarization(pol.x(),pol.y(),pol.z());
148
149 dp = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
150 dp->SetPolarization(-pol.x(),-pol.y(),-pol.z());
152
153 // Kill the incident positron
154 //
156 return &fParticleChange;
157}
158
159//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fAnnihilation
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetPolarization(G4double polX, G4double polY, G4double polZ)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void InitializeForPostStep(const G4Track &)
void AddSecondary(G4DynamicParticle *aParticle)
static G4Positron * Positron()
Definition: G4Positron.cc:94
Definition: G4Step.hh:78
G4double GetWeight() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetIntegral(G4bool val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetBuildTableFlag(G4bool val)
void SetEmModel(G4VEmModel *, G4int index=1)
G4VEmModel * EmModel(G4int index=1)
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4double MaxKinEnergy() const
G4double MinKinEnergy() const
void SetStartFromNullFlag(G4bool val)
G4ParticleChangeForGamma fParticleChange
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void SetNumberOfSecondaries(G4int totSecondaries)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
virtual void InitialiseProcess(const G4ParticleDefinition *)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
G4eplusAnnihilation(const G4String &name="annihil")
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)