Geant4 11.4.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#include "G4Log.hh"
69
70namespace
71{
72 constexpr G4double fTauPara = 0.12*CLHEP::ns;
73 constexpr G4double fTauOrto = 142.*CLHEP::ns;
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
94{
95 delete f2GammaAtRestModel;
96 delete f3GammaAtRestModel;
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118{
119 if (!isInitialised) {
120 isInitialised = true;
121 if(nullptr == EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
124 AddEmModel(1, EmModel(0));
125 }
126 auto param = G4EmParameters::Instance();
127
128 // AtRest models should be chosen only once
129 if (nullptr == f2GammaAtRestModel) {
130 auto type = param->PositronAtRestModelType();
131 if (type == fAllisonPositronium) {
132 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
133 } else if (type == fOrePowell) {
134 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
135 f3GammaAtRestModel = new G4OrePowellAtRestModel();
136 } else if (type == fOrePowellPolar) {
137 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
138 f3GammaAtRestModel = new G4PolarizedOrePowellAtRestModel();
139 } else {
140 f2GammaAtRestModel = new G4SimplePositronAtRestModel();
141 }
142 }
143 // Check that entanglement is switched on
144 // It may be set by the UI command "/process/em/QuantumEntanglement true".
145 fEntangled = param->QuantumEntanglement();
146 fApplyCuts = param->ApplyCuts();
147}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150
152{}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155
157 const G4Step& step)
158{
159 // positron at rest should be killed
160 fParticleChange.InitializeForPostStep(track);
161 fParticleChange.SetProposedKineticEnergy(0.);
162 fParticleChange.ProposeTrackStatus(fStopAndKill);
163
164 auto couple = step.GetPreStepPoint()->GetMaterialCutsCouple();
165 DefineMaterial(couple);
166
167 G4double gammaCut = GetGammaEnergyCut();
168
169 // apply cuts
170 if (fApplyCuts && gammaCut > CLHEP::electron_mass_c2) {
171 fParticleChange.ProposeLocalEnergyDeposit(2*CLHEP::electron_mass_c2);
172 return &fParticleChange;
173 }
174
175 // sample secondaries
176 secParticles.clear();
177 G4double edep = 0.0;
178 G4double ltime;
179 if (nullptr != f3GammaAtRestModel &&
180 G4UniformRand() < currentMaterial->GetIonisation()->GetOrtoPositroniumFraction()) {
181 f3GammaAtRestModel->SampleSecondaries(secParticles, edep, couple->GetMaterial());
182 ltime = fTauOrto;
183 } else {
184 f2GammaAtRestModel->SampleSecondaries(secParticles, edep, couple->GetMaterial());
185 ltime = fTauPara;
186 }
187
188 // define new weight for primary and secondaries
189 G4double weight = fParticleChange.GetParentWeight();
190 std::size_t num0 = secParticles.size();
191
192 // Russian roulette
193 if (nullptr != biasManager) {
194 G4int idx = couple->GetIndex();
195 if (biasManager->SecondaryBiasingRegion(idx) &&
196 !biasManager->GetDirectionalSplitting()) {
197 G4VEmModel* mod = EmModel(0);
198 G4double eloss = 0.0;
199 weight *= biasManager->ApplySecondaryBiasing(secParticles, track, mod,
200 &fParticleChange, eloss,
201 idx, gammaCut);
202 edep += eloss;
203 }
204 }
205
206 // save secondaries
207 std::size_t num = secParticles.size();
208
209 // Prepare a shared pointer only for two first gamma. If it is used, the
210 // shared pointer is copied into the tracks through G4EntanglementAuxInfo.
211 // This ensures the clip board lasts until both tracks are destroyed.
212 // It is assumed that 2 first secondary particles are the most energetic gamma
213 std::shared_ptr<G4eplusAnnihilationEntanglementClipBoard> clipBoard;
214 if (fEntangled && num >= 2) {
215 clipBoard = std::make_shared<G4eplusAnnihilationEntanglementClipBoard>();
216 clipBoard->SetParentParticleDefinition(track.GetDefinition());
217 }
218
219 if (num > 0) {
220 const G4double time = track.GetGlobalTime() - ltime*G4Log(G4UniformRand());
221 const G4ThreeVector& pos = track.GetPosition();
222 auto touch = track.GetTouchableHandle();
223 for (std::size_t i=0; i<num; ++i) {
225 G4Track* t = new G4Track(dp, time, pos);
226 t->SetTouchableHandle(touch);
227 if (fEntangled && i < 2) {
228 // entangledgammagamma is only true when there are only two gammas
229 // (See code above where entangledgammagamma is calculated.)
230 if (nullptr != clipBoard) {
231 if (i == 0) { // First gamma
232 clipBoard->SetTrackA(t);
233 } else if (i == 1) { // Second gamma
234 clipBoard->SetTrackB(t);
235 }
236 }
238 (fEntanglementModelID, new G4EntanglementAuxInfo(clipBoard));
239 }
240 if (nullptr != biasManager) {
241 t->SetWeight(weight * biasManager->GetWeight((G4int)i));
242 } else {
243 t->SetWeight(weight);
244 }
245 pParticleChange->AddSecondary(t);
246
247 // define type of secondary
248 if (i < num0) {
250 }
251 else {
253 }
254 }
255 }
256 fParticleChange.ProposeLocalEnergyDeposit(edep);
257 return &fParticleChange;
258}
259
260//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
261
262void G4eplusAnnihilation::ProcessDescription(std::ostream& out) const
263{
264 out << " Positron annihilation";
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
@ fAnnihilation
@ fEmDecreasing
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition G4Log.hh:169
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