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
G4AnnihiToMuPair.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// ------------ G4AnnihiToMuPair physics process ------
29// by H.Burkhardt, S. Kelner and R. Kokoulin, November 2002
30// -----------------------------------------------------------------------------
31//
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......//
33//
34// 27.01.03 : first implementation (hbu)
35// 04.02.03 : cosmetic simplifications (mma)
36// 25.10.04 : migrade to new interfaces of ParticleChange (vi)
37// 28.02.18 : cross section now including SSS threshold factor
38//
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40
41#include "G4AnnihiToMuPair.hh"
42
43#include "G4ios.hh"
44#include "Randomize.hh"
46#include "G4SystemOfUnits.hh"
47
48#include "G4Positron.hh"
49#include "G4MuonPlus.hh"
50#include "G4MuonMinus.hh"
51#include "G4TauPlus.hh"
52#include "G4TauMinus.hh"
53#include "G4Material.hh"
54#include "G4Step.hh"
55#include "G4LossTableManager.hh"
56#include "G4Exp.hh"
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60using namespace std;
61
63 G4ProcessType type):G4VDiscreteProcess (processName, type)
64{
65 //e+ Energy threshold
66 if(processName == "AnnihiToTauPair") {
68 part1 = G4TauPlus::TauPlus();
69 part2 = G4TauMinus::TauMinus();
70 fInfo = "e+e->tau+tau-";
71 } else {
73 part1 = G4MuonPlus::MuonPlus();
74 part2 = G4MuonMinus::MuonMinus();
75 }
76 fMass = part1->GetPDGMass();
77 fLowEnergyLimit =
78 2.*fMass*fMass/CLHEP::electron_mass_c2 - CLHEP::electron_mass_c2;
79
80 //model is ok up to 1000 TeV due to neglected Z-interference
81 fHighEnergyLimit = 1000.*TeV;
82
83 fCurrentSigma = 0.0;
84 fCrossSecFactor = 1.;
86 fManager->Register(this);
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
92{
93 fManager->DeRegister(this);
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
99{
100 return ( &particle == G4Positron::Positron() );
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104
106{
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111
113// Set the factor to artificially increase the cross section
114{
115 fCrossSecFactor = fac;
116 //G4cout << "The cross section for AnnihiToMuPair is artificially "
117 // << "increased by the CrossSecFactor=" << fCrossSecFactor << G4endl;
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
121
123// Calculates the microscopic cross section in GEANT4 internal units.
124// It gives a good description from threshold to 1000 GeV
125{
126 G4double rmuon = CLHEP::elm_coupling/fMass; //classical particle radius
127 G4double sig0 = CLHEP::pi*rmuon*rmuon/3.; //constant in crossSection
128 const G4double pial = CLHEP::pi*CLHEP::fine_structure_const; // pi * alphaQED
129
130 if (e <= fLowEnergyLimit) return 0.0;
131
132 const G4double xi = fLowEnergyLimit/e;
133 const G4double piaxi = pial * std::sqrt(xi);
134 G4double sigma = sig0 * xi * (1. + xi*0.5);
135 //G4cout << "### xi= " << xi << " piaxi=" << piaxi << G4endl;
136
137 // argument of the exponent below 0.1 or above 10
138 // Sigma per electron * number of electrons per atom
139 if(xi <= 1.0 - 100*piaxi*piaxi) {
140 sigma *= std::sqrt(1.0 - xi);
141 } else if( xi >= 1.0 - 0.01*piaxi*piaxi) {
142 sigma *= piaxi;
143 } else {
144 sigma *= piaxi/(1. - G4Exp( -piaxi/std::sqrt(1-xi) ));
145 }
146 //G4cout << "### sigma= " << sigma << G4endl;
147 return sigma;
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151
153 const G4double Z)
154{
155 return ComputeCrossSectionPerElectron(energy)*Z;
156}
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159
161 const G4Material* aMaterial)
162{
164}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167
170// returns the positron mean free path in GEANT4 internal units
171{
172 const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
173 G4double energy = aDynamicPositron->GetTotalEnergy();
174 const G4Material* aMaterial = aTrack.GetMaterial();
175
176 // cross section before step
177 fCurrentSigma = CrossSectionPerVolume(energy, aMaterial);
178
179 // increase the CrossSection by CrossSecFactor (default 1)
180 return (fCurrentSigma > 0.0) ? 1.0/(fCurrentSigma*fCrossSecFactor) : DBL_MAX;
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184
186 const G4Step& aStep)
187//
188// generation of e+e- -> mu+mu-
189//
190{
192
193 // current Positron energy and direction, return if energy too low
194 const G4DynamicParticle *aDynamicPositron = aTrack.GetDynamicParticle();
195 const G4double Mele = CLHEP::electron_mass_c2;
196 G4double Epos = aDynamicPositron->GetTotalEnergy();
197 G4double xs = CrossSectionPerVolume(Epos, aTrack.GetMaterial());
198
199 // test of cross section
200 if(xs > 0.0 && fCurrentSigma*G4UniformRand() > xs)
201 {
202 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
203 }
204
205 const G4ThreeVector PosiDirection = aDynamicPositron->GetMomentumDirection();
206 G4double xi = fLowEnergyLimit/Epos; // xi is always less than 1,
207 // goes to 0 at high Epos
208
209 // generate cost; probability function 1+cost**2 at high Epos
210 //
211 G4double cost;
212 do { cost = 2.*G4UniformRand()-1.; }
213 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
214 while (2.*G4UniformRand() > 1.+xi+cost*cost*(1.-xi) );
215 G4double sint = sqrt(1.-cost*cost);
216
217 // generate phi
218 //
219 G4double phi = 2.*CLHEP::pi*G4UniformRand();
220
221 G4double Ecm = std::sqrt(0.5*Mele*(Epos+Mele));
222 G4double Pcm = std::sqrt(Ecm*Ecm - fMass*fMass);
223 G4double beta = std::sqrt((Epos-Mele)/(Epos+Mele));
224 G4double gamma = Ecm/Mele;
225 G4double Pt = Pcm*sint;
226
227 // energy and momentum of the muons in the Lab
228 //
229 G4double EmuPlus = gamma*(Ecm + cost*beta*Pcm);
230 G4double EmuMinus = gamma*(Ecm - cost*beta*Pcm);
231 G4double PmuPlusZ = gamma*(beta*Ecm + cost*Pcm);
232 G4double PmuMinusZ = gamma*(beta*Ecm - cost*Pcm);
233 G4double PmuPlusX = Pt*std::cos(phi);
234 G4double PmuPlusY = Pt*std::sin(phi);
235 G4double PmuMinusX =-PmuPlusX;
236 G4double PmuMinusY =-PmuPlusY;
237 // absolute momenta
238 G4double PmuPlus = std::sqrt(Pt*Pt+PmuPlusZ *PmuPlusZ );
239 G4double PmuMinus = std::sqrt(Pt*Pt+PmuMinusZ*PmuMinusZ);
240
241 // mu+ mu- directions for Positron in z-direction
242 //
244 MuPlusDirection(PmuPlusX/PmuPlus, PmuPlusY/PmuPlus, PmuPlusZ/PmuPlus);
246 MuMinusDirection(PmuMinusX/PmuMinus,PmuMinusY/PmuMinus,PmuMinusZ/PmuMinus);
247
248 // rotate to actual Positron direction
249 //
250 MuPlusDirection.rotateUz(PosiDirection);
251 MuMinusDirection.rotateUz(PosiDirection);
252
254
255 // create G4DynamicParticle object for the particle1
256 G4DynamicParticle* aParticle1 =
257 new G4DynamicParticle(part1, MuPlusDirection, EmuPlus-fMass);
258 aParticleChange.AddSecondary(aParticle1);
259 // create G4DynamicParticle object for the particle2
260 G4DynamicParticle* aParticle2 =
261 new G4DynamicParticle(part2, MuMinusDirection, EmuMinus-fMass);
262 aParticleChange.AddSecondary(aParticle2);
263
264 // Kill the incident positron
265 //
268
269 return &aParticleChange;
270}
271
272//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
273
275{
276 G4String comments = fInfo + " annihilation, atomic e- at rest, SubType=";
277 G4cout << G4endl << GetProcessName() << ": " << comments
279 G4cout << " threshold at " << fLowEnergyLimit/CLHEP::GeV << " GeV"
280 << " good description up to "
281 << fHighEnergyLimit/CLHEP::TeV << " TeV for all Z." << G4endl;
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fAnnihilationToTauTau
@ fAnnihilationToMuMu
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4ForceCondition
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *) override
G4double ComputeCrossSectionPerElectron(const G4double energy)
~G4AnnihiToMuPair() override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4bool IsApplicable(const G4ParticleDefinition &) override
void SetCrossSecFactor(G4double fac)
G4double ComputeCrossSectionPerAtom(const G4double energy, const G4double Z)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4AnnihiToMuPair(const G4String &processName="AnnihiToMuPair", G4ProcessType type=fElectromagnetic)
G4double CrossSectionPerVolume(G4double positronEnergy, const G4Material *)
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalEnergy() const
static G4LossTableManager * Instance()
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4double GetTotNbOfElectPerVolume() const
Definition: G4Material.hh:207
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
static G4Positron * Positron()
Definition: G4Positron.cc:93
Definition: G4Step.hh:62
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:134
static G4TauPlus * TauPlus()
Definition: G4TauPlus.cc:133
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define DBL_MAX
Definition: templates.hh:62