source: trunk/source/processes/electromagnetic/highenergy/src/G4AnnihiToMuPair.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 9.7 KB
RevLine 
[819]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//
[961]27// $Id: G4AnnihiToMuPair.cc,v 1.5 2008/10/16 14:29:48 vnivanch Exp $
[1007]28// GEANT4 tag $Name: geant4-09-02 $
[819]29//
30// ------------ G4AnnihiToMuPair physics process ------
31// by H.Burkhardt, S. Kelner and R. Kokoulin, November 2002
32// -----------------------------------------------------------------------------
33//
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......//
35//
36// 27.01.03 : first implementation (hbu)
37// 04.02.03 : cosmetic simplifications (mma)
38// 25.10.04 : migrade to new interfaces of ParticleChange (vi)
39//
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
42#include "G4AnnihiToMuPair.hh"
43
44#include "G4ios.hh"
45#include "Randomize.hh"
46
47#include "G4Positron.hh"
48#include "G4MuonPlus.hh"
49#include "G4MuonMinus.hh"
50#include "G4Material.hh"
51#include "G4Step.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54
55using namespace std;
56
57G4AnnihiToMuPair::G4AnnihiToMuPair(const G4String& processName,
58 G4ProcessType type):G4VDiscreteProcess (processName, type)
59{
60 //e+ Energy threshold
61 const G4double Mu_massc2 = G4MuonPlus::MuonPlus()->GetPDGMass();
62 LowestEnergyLimit = 2*Mu_massc2*Mu_massc2/electron_mass_c2 - electron_mass_c2;
63
64 //modele ok up to 1000 TeV due to neglected Z-interference
65 HighestEnergyLimit = 1000*TeV;
66
67 CrossSecFactor = 1.;
[961]68 SetProcessSubType(6);
69
[819]70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74G4AnnihiToMuPair::~G4AnnihiToMuPair() // (empty) destructor
75{ }
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
79G4bool G4AnnihiToMuPair::IsApplicable(const G4ParticleDefinition& particle)
80{
81 return ( &particle == G4Positron::Positron() );
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
86void G4AnnihiToMuPair::BuildPhysicsTable(const G4ParticleDefinition&)
87// Build cross section and mean free path tables
88//here no tables, just calling PrintInfoDefinition
89{
90 PrintInfoDefinition();
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94
95void G4AnnihiToMuPair::SetCrossSecFactor(G4double fac)
96// Set the factor to artificially increase the cross section
97{
98 CrossSecFactor = fac;
99 G4cout << "The cross section for AnnihiToMuPair is artificially "
100 << "increased by the CrossSecFactor=" << CrossSecFactor << G4endl;
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104
105G4double G4AnnihiToMuPair::ComputeCrossSectionPerAtom(G4double Epos, G4double Z)
106// Calculates the microscopic cross section in GEANT4 internal units.
107// It gives a good description from threshold to 1000 GeV
108{
109 static const G4double Mmuon = G4MuonPlus::MuonPlus()->GetPDGMass();
110 static const G4double Rmuon = elm_coupling/Mmuon; //classical particle radius
111 static const G4double Sig0 = pi*Rmuon*Rmuon/3.; //constant in crossSection
112
113 G4double CrossSection = 0.;
114 if (Epos < LowestEnergyLimit) return CrossSection;
115
116 G4double xi = LowestEnergyLimit/Epos;
117 G4double SigmaEl = Sig0*xi*(1.+xi/2.)*sqrt(1.-xi); // per electron
118 CrossSection = SigmaEl*Z; // number of electrons per atom
119 CrossSection *= CrossSecFactor; //increase the CrossSection by (default 1)
120 return CrossSection;
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
125G4double G4AnnihiToMuPair::GetMeanFreePath(const G4Track& aTrack,
126 G4double, G4ForceCondition*)
127
128// returns the positron mean free path in GEANT4 internal units
129
130{
131 const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
132 G4double PositronEnergy = aDynamicPositron->GetKineticEnergy()
133 +electron_mass_c2;
134 G4Material* aMaterial = aTrack.GetMaterial();
135 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
136 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
137
138 G4double SIGMA = 0 ;
139
140 for ( size_t i=0 ; i < aMaterial->GetNumberOfElements() ; i++ )
141 {
142 G4double AtomicZ = (*theElementVector)[i]->GetZ();
143 SIGMA += NbOfAtomsPerVolume[i] *
144 ComputeCrossSectionPerAtom(PositronEnergy,AtomicZ);
145 }
146 return SIGMA > DBL_MIN ? 1./SIGMA : DBL_MAX;
147}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150
151G4VParticleChange* G4AnnihiToMuPair::PostStepDoIt(const G4Track& aTrack,
152 const G4Step& aStep)
153//
154// generation of e+e- -> mu+mu-
155//
156{
157
158 aParticleChange.Initialize(aTrack);
159 static const G4double Mele=electron_mass_c2;
160 static const G4double Mmuon=G4MuonPlus::MuonPlus()->GetPDGMass();
161
162 // current Positron energy and direction, return if energy too low
163 const G4DynamicParticle *aDynamicPositron = aTrack.GetDynamicParticle();
164 G4double Epos = aDynamicPositron->GetKineticEnergy()+Mele;
165
166 if (Epos < LowestEnergyLimit)
167 { G4cout
168 << "error in G4AnnihiToMuPair::PostStepDoIt called with energy below"
169 " threshold Epos= "
170 << Epos << G4endl; // shoud never happen
171 G4Exception(10);
172 }
173
174 if (Epos < LowestEnergyLimit)
175 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
176
177 G4ParticleMomentum PositronDirection =
178 aDynamicPositron->GetMomentumDirection();
179 G4double xi = LowestEnergyLimit/Epos; // xi is always less than 1,
180 // goes to 0 at high Epos
181
182 // generate cost
183 //
184 G4double cost;
185 do cost = 2.*G4UniformRand()-1.;
186 while (2.*G4UniformRand() > 1.+xi+cost*cost*(1.-xi) );
187 //1+cost**2 at high Epos
188 G4double sint = sqrt(1.-cost*cost);
189
190 // generate phi
191 //
192 G4double phi=2.*pi*G4UniformRand();
193
194 G4double Ecm = sqrt(0.5*Mele*(Epos+Mele));
195 G4double Pcm = sqrt(Ecm*Ecm-Mmuon*Mmuon);
196 G4double beta = sqrt((Epos-Mele)/(Epos+Mele));
197 G4double gamma = Ecm/Mele; // =sqrt((Epos+Mele)/(2.*Mele));
198 G4double Pt = Pcm*sint;
199
200 // energy and momentum of the muons in the Lab
201 //
202 G4double EmuPlus = gamma*( Ecm+cost*beta*Pcm);
203 G4double EmuMinus = gamma*( Ecm-cost*beta*Pcm);
204 G4double PmuPlusZ = gamma*(beta*Ecm+cost* Pcm);
205 G4double PmuMinusZ = gamma*(beta*Ecm-cost* Pcm);
206 G4double PmuPlusX = Pt*cos(phi);
207 G4double PmuPlusY = Pt*sin(phi);
208 G4double PmuMinusX =-Pt*cos(phi);
209 G4double PmuMinusY =-Pt*sin(phi);
210 // absolute momenta
211 G4double PmuPlus = sqrt(Pt*Pt+PmuPlusZ *PmuPlusZ );
212 G4double PmuMinus = sqrt(Pt*Pt+PmuMinusZ*PmuMinusZ);
213
214 // mu+ mu- directions for Positron in z-direction
215 //
216 G4ThreeVector
217 MuPlusDirection ( PmuPlusX/PmuPlus, PmuPlusY/PmuPlus, PmuPlusZ/PmuPlus );
218 G4ThreeVector
219 MuMinusDirection(PmuMinusX/PmuMinus,PmuMinusY/PmuMinus,PmuMinusZ/PmuMinus);
220
221 // rotate to actual Positron direction
222 //
223 MuPlusDirection.rotateUz(PositronDirection);
224 MuMinusDirection.rotateUz(PositronDirection);
225
226 aParticleChange.SetNumberOfSecondaries(2);
227 // create G4DynamicParticle object for the particle1
228 G4DynamicParticle* aParticle1= new G4DynamicParticle(
229 G4MuonPlus::MuonPlus(),MuPlusDirection,EmuPlus-Mmuon);
230 aParticleChange.AddSecondary(aParticle1);
231 // create G4DynamicParticle object for the particle2
232 G4DynamicParticle* aParticle2= new G4DynamicParticle(
233 G4MuonMinus::MuonMinus(),MuMinusDirection,EmuMinus-Mmuon);
234 aParticleChange.AddSecondary(aParticle2);
235
236 // Kill the incident positron
237 //
238 aParticleChange.ProposeEnergy(0.);
239 aParticleChange.ProposeTrackStatus(fStopAndKill);
240
241 return &aParticleChange;
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245
246void G4AnnihiToMuPair::PrintInfoDefinition()
247{
[961]248 G4String comments ="e+e->mu+mu- annihilation, atomic e- at rest, SubType=.";
249 G4cout << G4endl << GetProcessName() << ": " << comments
250 << GetProcessSubType() << G4endl;
251 G4cout << " threshold at " << LowestEnergyLimit/GeV << " GeV"
[819]252 << " good description up to "
253 << HighestEnergyLimit/TeV << " TeV for all Z." << G4endl;
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.