Geant4-11
All Data Structures 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
70{
74 SetBuildTableFlag(false);
78 enableAtRestDoIt = true;
80 fEntanglementModelID = G4PhysicsModelCatalog::GetModelID("model_GammaGammaEntanglement");
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86{}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
91{
92 return (&p == G4Positron::Positron());
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
99{
101 return 0.0;
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
107{
108 if(!isInitialised) {
109 isInitialised = true;
110 if(nullptr == EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
113 AddEmModel(1, EmModel(0));
114 }
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118
120{}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123
125 const G4Step& step)
126// Performs the e+ e- annihilation when both particles are assumed at rest.
127{
129
131 size_t idx = CurrentMaterialCutsCoupleIndex();
132 G4double ene(0.0);
133 G4VEmModel* model = SelectModel(ene, idx);
134
135 // define new weight for primary and secondaries
137
138 // sample secondaries
139 secParticles.clear();
140 G4double gammaCut = GetGammaEnergyCut();
142 track.GetDynamicParticle(), gammaCut);
143
144 G4int num0 = secParticles.size();
145
146 // splitting or Russian roulette
147 if(biasManager) {
149 G4double eloss = 0.0;
151 secParticles, track, model, &fParticleChange, eloss,
152 idx, gammaCut, step.GetPostStepPoint()->GetSafety());
153 if(eloss > 0.0) {
156 }
157 }
158 }
159
160 // save secondaries
161 G4int num = secParticles.size();
162
163 // Check that entanglement is switched on... (the following flag is
164 // set by /process/em/QuantumEntanglement).
166 // ...and that we have two gammas with both gammas' energies above
167 // gammaCut (entanglement is only programmed for e+ e- -> gamma gamma).
168 G4bool entangledgammagamma = false;
169 if (entangled) {
170 if (num == 2) {
171 entangledgammagamma = true;
172 for (const auto* p: secParticles) {
173 if (p->GetDefinition() != theGamma ||
174 p->GetKineticEnergy() < gammaCut) {
175 entangledgammagamma = false;
176 }
177 }
178 }
179 }
180
181 // Prepare a shared pointer for psossible use below. If it is used, the
182 // shared pointer is copied into the tracks through G4EntanglementAuxInfo.
183 // This ensures the clip board lasts until both tracks are destroyed.
184 std::shared_ptr<G4eplusAnnihilationEntanglementClipBoard> clipBoard;
185 if (entangledgammagamma) {
186 clipBoard = std::make_shared<G4eplusAnnihilationEntanglementClipBoard>();
187 clipBoard->SetParentParticleDefinition(track.GetDefinition());
188 }
189
190 if(num > 0) {
191
194 G4double time = track.GetGlobalTime();
195
196 for (G4int i=0; i<num; ++i) {
197 if (secParticles[i]) {
200 G4double e = dp->GetKineticEnergy();
201 G4bool good = true;
202 if(ApplyCuts()) {
203 if (p == theGamma) {
204 if (e < gammaCut) { good = false; }
205 } else if (p == theElectron) {
206 if (e < GetElectronEnergyCut()) { good = false; }
207 }
208 // added secondary if it is good
209 }
210 if (good) {
211 G4Track* t = new G4Track(dp, time, track.GetPosition());
213 if (entangledgammagamma) {
214 // entangledgammagamma is only true when there are only two gammas
215 // (See code above where entangledgammagamma is calculated.)
216 if (i == 0) { // First gamma
217 clipBoard->SetTrackA(t);
218 } else if (i == 1) { // Second gamma
219 clipBoard->SetTrackB(t);
220 }
223 }
224 if (biasManager) {
225 t->SetWeight(weight * biasManager->GetWeight(i));
226 } else {
227 t->SetWeight(weight);
228 }
230
231 // define type of secondary
233 else if(i < num0) {
234 if(p == theGamma) {
236 } else {
238 }
239 } else {
241 }
242 /*
243 G4cout << "Secondary(post step) has weight " << t->GetWeight()
244 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
245 << GetProcessName() << " fluoID= " << fluoID
246 << " augerID= " << augerID <<G4endl;
247 */
248 } else {
249 delete dp;
250 edep += e;
251 }
252 }
253 }
255 }
256 return &fParticleChange;
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260
261void G4eplusAnnihilation::ProcessDescription(std::ostream& out) const
262{
263 out << " Positron annihilation";
265}
266
267//....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:205
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:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
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:359
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4ParticleDefinition * theElectron
void ProcessDescription(std::ostream &) const override
const G4ParticleDefinition * theGamma
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
const char * name(G4int ptype)