Geant4 10.7.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4eplusAnnihilation
33//
34// Author: Vladimir Ivanchenko on base of Michel Maire code
35//
36// Creation date: 02.08.2004
37//
38// Modified by Michel Maire, Vladimir Ivanchenko and Daren Sawkey
39
40//
41// -------------------------------------------------------------------
42//
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45
49#include "G4Gamma.hh"
50#include "G4Positron.hh"
52#include "G4EmBiasingManager.hh"
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
56using namespace std;
57
59 : G4VEmProcess(name), isInitialised(false)
60{
61 theGamma = G4Gamma::Gamma();
62 SetIntegral(true);
63 SetBuildTableFlag(false);
65 SetSecondaryParticle(theGamma);
67 enableAtRestDoIt = true;
69}
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
74{}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
79{
80 return (&p == G4Positron::Positron());
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
87{
89 return 0.0;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95{
96 if(!isInitialised) {
97 isInitialised = true;
98 if(!EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
101 AddEmModel(1, EmModel(0));
102 }
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
108{}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111
113 const G4Step& step)
114// Performs the e+ e- annihilation when both particles are assumed at rest.
115{
117
119 size_t idx = CurrentMaterialCutsCoupleIndex();
120 G4double ene(0.0);
121 G4VEmModel* model = SelectModel(ene, idx);
122
123 // define new weight for primary and secondaries
125
126 // sample secondaries
127 secParticles.clear();
128 G4double gammaCut = GetGammaEnergyCut();
130 track.GetDynamicParticle(), gammaCut);
131
132 G4int num0 = secParticles.size();
133
134 // splitting or Russian roulette
135 if(biasManager) {
137 G4double eloss = 0.0;
139 secParticles, track, model, &fParticleChange, eloss,
140 idx, gammaCut, step.GetPostStepPoint()->GetSafety());
141 if(eloss > 0.0) {
144 }
145 }
146 }
147 // save secondaries
148 G4int num = secParticles.size();
149 if(num > 0) {
150
153 G4double time = track.GetGlobalTime();
154
155 for (G4int i=0; i<num; ++i) {
156 if (secParticles[i]) {
159 G4double e = dp->GetKineticEnergy();
160 G4bool good = true;
161 if(ApplyCuts()) {
162 if (p == theGamma) {
163 if (e < gammaCut) { good = false; }
164 } else if (p == theElectron) {
165 if (e < GetElectronEnergyCut()) { good = false; }
166 }
167 // added secondary if it is good
168 }
169 if (good) {
170 G4Track* t = new G4Track(dp, time, track.GetPosition());
172 if (biasManager) {
173 t->SetWeight(weight * biasManager->GetWeight(i));
174 } else {
175 t->SetWeight(weight);
176 }
178
179 // define type of secondary
181 else if(i < num0) {
182 if(p == theGamma) {
184 } else {
186 }
187 } else {
189 }
190 /*
191 G4cout << "Secondary(post step) has weight " << t->GetWeight()
192 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
193 << GetProcessName() << " fluoID= " << fluoID
194 << " augerID= " << augerID <<G4endl;
195 */
196 } else {
197 delete dp;
198 edep += e;
199 }
200 }
201 }
203 }
204 return &fParticleChange;
205}
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208
209void G4eplusAnnihilation::ProcessDescription(std::ostream& out) const
210{
211 out << " Positron annihilation";
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fAnnihilation
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
G4double GetWeight(G4int i)
G4bool SecondaryBiasingRegion(G4int coupleIdx)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void InitializeForPostStep(const G4Track &)
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4double GetSafety() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetCreatorModelIndex(G4int idx)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4bool ApplyCuts() const
G4double GetGammaEnergyCut()
void SetIntegral(G4bool val)
G4int mainSecondaries
G4EmBiasingManager * biasManager
G4VEmModel * EmModel(size_t index=0) const
void SetBuildTableFlag(G4bool val)
G4double GetElectronEnergyCut()
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetSecondaryParticle(const G4ParticleDefinition *p)
virtual void ProcessDescription(std::ostream &outFile) const override
G4double MaxKinEnergy() const
G4VEmModel * SelectModel(G4double kinEnergy, size_t index)
G4double MinKinEnergy() const
std::vector< G4DynamicParticle * > secParticles
const G4ParticleDefinition * theElectron
void SetStartFromNullFlag(G4bool val)
const G4MaterialCutsCouple * MaterialCutsCouple() const
G4ParticleChangeForGamma fParticleChange
size_t CurrentMaterialCutsCoupleIndex() const
G4double GetParentWeight() const
G4double GetLocalEnergyDeposit() const
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
virtual void ProcessDescription(std::ostream &) const override
virtual void InitialiseProcess(const G4ParticleDefinition *) override
G4eplusAnnihilation(const G4String &name="annihil")
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition) override
virtual G4bool IsApplicable(const G4ParticleDefinition &p) final
virtual void StreamProcessInfo(std::ostream &outFile) const override
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData) override