Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "G4Electron.hh"
58#include "G4Positron.hh"
60#include "G4EmBiasingManager.hh"
63#include "G4EmParameters.hh"
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
69 : G4VEmProcess(name)
70{
71 theGamma = G4Gamma::Gamma();
72 theElectron = G4Electron::Electron();
74 SetBuildTableFlag(false);
76 SetSecondaryParticle(theGamma);
78 enableAtRestDoIt = true;
80 fEntanglementModelID = G4PhysicsModelCatalog::GetModelID("model_GammaGammaEntanglement");
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90{
91 return (&p == G4Positron::Positron());
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95
98{
100 return 0.0;
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104
106{
107 if(!isInitialised) {
108 isInitialised = true;
109 if(nullptr == EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
112 AddEmModel(1, EmModel(0));
113 }
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117
119{}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122
124 const G4Step& step)
125// Performs the e+ e- annihilation when both particles are assumed at rest.
126{
128
131 G4double ene(0.0);
132 G4VEmModel* model = SelectModel(ene, idx);
133
134 // define new weight for primary and secondaries
136
137 // sample secondaries
138 secParticles.clear();
139 G4double gammaCut = GetGammaEnergyCut();
141 track.GetDynamicParticle(), gammaCut);
142
143 G4int num0 = (G4int)secParticles.size();
144
145 // splitting or Russian roulette
146 if(biasManager) {
148 G4double eloss = 0.0;
150 secParticles, track, model, &fParticleChange, eloss,
151 idx, gammaCut, step.GetPostStepPoint()->GetSafety());
152 if(eloss > 0.0) {
155 }
156 }
157 }
158
159 // save secondaries
160 G4int num = (G4int)secParticles.size();
161
162 // Check that entanglement is switched on... (the following flag is
163 // set by /process/em/QuantumEntanglement).
165 // ...and that we have two gammas with both gammas' energies above
166 // gammaCut (entanglement is only programmed for e+ e- -> gamma gamma).
167 G4bool entangledgammagamma = false;
168 if (entangled) {
169 if (num == 2) {
170 entangledgammagamma = true;
171 for (const auto* p: secParticles) {
172 if (p->GetDefinition() != theGamma ||
173 p->GetKineticEnergy() < gammaCut) {
174 entangledgammagamma = false;
175 }
176 }
177 }
178 }
179
180 // Prepare a shared pointer for psossible use below. If it is used, the
181 // shared pointer is copied into the tracks through G4EntanglementAuxInfo.
182 // This ensures the clip board lasts until both tracks are destroyed.
183 std::shared_ptr<G4eplusAnnihilationEntanglementClipBoard> clipBoard;
184 if (entangledgammagamma) {
185 clipBoard = std::make_shared<G4eplusAnnihilationEntanglementClipBoard>();
186 clipBoard->SetParentParticleDefinition(track.GetDefinition());
187 }
188
189 if(num > 0) {
190
193 G4double time = track.GetGlobalTime();
194
195 for (G4int i=0; i<num; ++i) {
196 if (secParticles[i]) {
199 G4double e = dp->GetKineticEnergy();
200 G4bool good = true;
201 if(ApplyCuts()) {
202 if (p == theGamma) {
203 if (e < gammaCut) { good = false; }
204 } else if (p == theElectron) {
205 if (e < GetElectronEnergyCut()) { good = false; }
206 }
207 // added secondary if it is good
208 }
209 if (good) {
210 G4Track* t = new G4Track(dp, time, track.GetPosition());
212 if (entangledgammagamma) {
213 // entangledgammagamma is only true when there are only two gammas
214 // (See code above where entangledgammagamma is calculated.)
215 if (i == 0) { // First gamma
216 clipBoard->SetTrackA(t);
217 } else if (i == 1) { // Second gamma
218 clipBoard->SetTrackB(t);
219 }
221 (fEntanglementModelID,new G4EntanglementAuxInfo(clipBoard));
222 }
223 if (biasManager) {
224 t->SetWeight(weight * biasManager->GetWeight(i));
225 } else {
226 t->SetWeight(weight);
227 }
229
230 // define type of secondary
232 else if(i < num0) {
233 if(p == theGamma) {
235 } else {
237 }
238 } else {
240 }
241 /*
242 G4cout << "Secondary(post step) has weight " << t->GetWeight()
243 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
244 << GetProcessName() << " fluoID= " << fluoID
245 << " augerID= " << augerID <<G4endl;
246 */
247 } else {
248 delete dp;
249 edep += e;
250 }
251 }
252 }
254 }
255 return &fParticleChange;
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
259
260void G4eplusAnnihilation::ProcessDescription(std::ostream& out) const
261{
262 out << " Positron annihilation";
264}
265
266//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fAnnihilation
@ fEmDecreasing
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
static G4Electron * Electron()
Definition: G4Electron.cc:93
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 G4EmParameters * Instance()
G4bool QuantumEntanglement() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void InitializeForPostStep(const G4Track &)
static G4int GetModelID(const G4int modelIndex)
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4double GetSafety() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:209
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
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()
G4int mainSecondaries
G4VEmModel * SelectModel(G4double kinEnergy, size_t)
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)
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)
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:363
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
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
~G4eplusAnnihilation() override