Geant4 11.3.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// Introduced Quantum Entanglement April 2021 John Allison
41// This is activated by /process/em/QuantumEntanglement
42// For e+e- -> gamma gamma, the gammas are "tagged" here
43// and must be "analysed" in a Compton scattering process - see, for
44// example, G4LivermorePolarizedComptonModel. Otherwise entanglement
45// has no effect even if activated.
46
47//
48// -------------------------------------------------------------------
49//
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52
56#include "G4Gamma.hh"
57#include "G4Positron.hh"
59#include "G4EmBiasingManager.hh"
66#include "G4EmParameters.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87{
88 delete f2GammaAtRestModel;
89 delete f3GammaAtRestModel;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109
111{
112 if (!isInitialised) {
113 isInitialised = true;
114 if(nullptr == EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
117 AddEmModel(1, EmModel(0));
118 }
119 auto param = G4EmParameters::Instance();
120
121 // AtRest models should be chosen only once
122 if (nullptr == f2GammaAtRestModel) {
123 auto type = param->PositronAtRestModelType();
124 if (type == fAllisonPositronium) {
125 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
126 } else if (type == fOrePowell) {
127 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
128 f3GammaAtRestModel = new G4OrePowellAtRestModel();
129 } else if (type == fOrePowellPolar) {
130 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
131 f3GammaAtRestModel = new G4PolarizedOrePowellAtRestModel();
132 } else {
133 f2GammaAtRestModel = new G4SimplePositronAtRestModel();
134 }
135 }
136 // Check that entanglement is switched on
137 // It may be set by the UI command "/process/em/QuantumEntanglement true".
138 fEntangled = param->QuantumEntanglement();
139 fApplyCuts = param->ApplyCuts();
140}
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
143
145{}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148
150 const G4Step& step)
151{
152 // positron at rest should be killed
153 fParticleChange.InitializeForPostStep(track);
154 fParticleChange.SetProposedKineticEnergy(0.);
155 fParticleChange.ProposeTrackStatus(fStopAndKill);
156
157 auto couple = step.GetPreStepPoint()->GetMaterialCutsCouple();
158 DefineMaterial(couple);
159
160 G4double gammaCut = GetGammaEnergyCut();
161
162 // apply cuts
163 if (fApplyCuts && gammaCut > CLHEP::electron_mass_c2) {
164 fParticleChange.ProposeLocalEnergyDeposit(2*CLHEP::electron_mass_c2);
165 return &fParticleChange;
166 }
167
168 // sample secondaries
169 secParticles.clear();
170 G4double edep = 0.0;
171 if (nullptr != f3GammaAtRestModel &&
172 G4UniformRand() < currentMaterial->GetIonisation()->GetOrtoPositroniumFraction()) {
173 f3GammaAtRestModel->SampleSecondaries(secParticles, edep, couple->GetMaterial());
174 } else {
175 f2GammaAtRestModel->SampleSecondaries(secParticles, edep, couple->GetMaterial());
176 }
177
178 // define new weight for primary and secondaries
179 G4double weight = fParticleChange.GetParentWeight();
180 std::size_t num0 = secParticles.size();
181
182 // Russian roulette
183 if (nullptr != biasManager) {
184 G4int idx = couple->GetIndex();
185 if (biasManager->SecondaryBiasingRegion(idx) &&
186 !biasManager->GetDirectionalSplitting()) {
187 G4VEmModel* mod = EmModel(0);
188 G4double eloss = 0.0;
189 weight *= biasManager->ApplySecondaryBiasing(secParticles, track, mod,
190 &fParticleChange, eloss,
191 idx, gammaCut);
192 edep += eloss;
193 }
194 }
195
196 // save secondaries
197 std::size_t num = secParticles.size();
198
199 // Prepare a shared pointer only for two first gamma. If it is used, the
200 // shared pointer is copied into the tracks through G4EntanglementAuxInfo.
201 // This ensures the clip board lasts until both tracks are destroyed.
202 // It is assumed that 2 first secondary particles are the most energetic gamma
203 std::shared_ptr<G4eplusAnnihilationEntanglementClipBoard> clipBoard;
204 if (fEntangled && num >= 2) {
205 clipBoard = std::make_shared<G4eplusAnnihilationEntanglementClipBoard>();
206 clipBoard->SetParentParticleDefinition(track.GetDefinition());
207 }
208
209 if (num > 0) {
210 const G4double time = track.GetGlobalTime();
211 const G4ThreeVector& pos = track.GetPosition();
212 auto touch = track.GetTouchableHandle();
213 for (std::size_t i=0; i<num; ++i) {
215 G4Track* t = new G4Track(dp, time, pos);
216 t->SetTouchableHandle(touch);
217 if (fEntangled && i < 2) {
218 // entangledgammagamma is only true when there are only two gammas
219 // (See code above where entangledgammagamma is calculated.)
220 if (nullptr != clipBoard) {
221 if (i == 0) { // First gamma
222 clipBoard->SetTrackA(t);
223 } else if (i == 1) { // Second gamma
224 clipBoard->SetTrackB(t);
225 }
226 }
228 (fEntanglementModelID, new G4EntanglementAuxInfo(clipBoard));
229 }
230 if (nullptr != biasManager) {
231 t->SetWeight(weight * biasManager->GetWeight((G4int)i));
232 } else {
233 t->SetWeight(weight);
234 }
235 pParticleChange->AddSecondary(t);
236
237 // define type of secondary
238 if (i < num0) {
240 }
241 else {
243 }
244 }
245 }
246 fParticleChange.ProposeLocalEnergyDeposit(edep);
247 return &fParticleChange;
248}
249
250//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
251
252void G4eplusAnnihilation::ProcessDescription(std::ostream& out) const
253{
254 out << " Positron annihilation";
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
@ fAnnihilation
@ fEmDecreasing
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
static G4EmParameters * Instance()
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4int GetModelID(const G4int modelIndex)
static G4Positron * Positron()
Definition G4Positron.cc:90
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4StepPoint * GetPreStepPoint() const
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition G4Track.cc:215
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4double GetGammaEnergyCut()
G4int mainSecondaries
G4VEmModel * EmModel(std::size_t index=0) const
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
G4EmBiasingManager * biasManager
void SetBuildTableFlag(G4bool val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetCrossSectionType(G4CrossSectionType val)
void ProcessDescription(std::ostream &outFile) const override
G4double MaxKinEnergy() const
G4double MinKinEnergy() const
std::vector< G4DynamicParticle * > secParticles
void SetStartFromNullFlag(G4bool val)
G4ParticleChangeForGamma fParticleChange
const G4Material * currentMaterial
G4bool enableAtRestDoIt
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
void ProcessDescription(std::ostream &) const override
void InitialiseProcess(const G4ParticleDefinition *) override
G4eplusAnnihilation(const G4String &name="annihil")
G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition) override
G4bool IsApplicable(const G4ParticleDefinition &p) final
void StreamProcessInfo(std::ostream &outFile) const override
G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData) override