Geant4  10.02.p02
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 // $Id: G4AnnihiToMuPair.cc 91869 2015-08-07 15:21:02Z gcosmo $
28 //
29 // ------------ G4AnnihiToMuPair physics process ------
30 // by H.Burkhardt, S. Kelner and R. Kokoulin, November 2002
31 // -----------------------------------------------------------------------------
32 //
33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......//
34 //
35 // 27.01.03 : first implementation (hbu)
36 // 04.02.03 : cosmetic simplifications (mma)
37 // 25.10.04 : migrade to new interfaces of ParticleChange (vi)
38 //
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40 
41 #include "G4AnnihiToMuPair.hh"
42 
43 #include "G4ios.hh"
44 #include "Randomize.hh"
45 #include "G4PhysicalConstants.hh"
46 #include "G4SystemOfUnits.hh"
47 
48 #include "G4Positron.hh"
49 #include "G4MuonPlus.hh"
50 #include "G4MuonMinus.hh"
51 #include "G4Material.hh"
52 #include "G4Step.hh"
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
56 using namespace std;
57 
59  G4ProcessType type):G4VDiscreteProcess (processName, type)
60 {
61  //e+ Energy threshold
62  const G4double Mu_massc2 = G4MuonPlus::MuonPlus()->GetPDGMass();
63  LowestEnergyLimit = 2*Mu_massc2*Mu_massc2/electron_mass_c2 - electron_mass_c2;
64 
65  //modele ok up to 1000 TeV due to neglected Z-interference
66  HighestEnergyLimit = 1000*TeV;
67 
68  CurrentSigma = 0.0;
69  CrossSecFactor = 1.;
71 
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75 
76 G4AnnihiToMuPair::~G4AnnihiToMuPair() // (empty) destructor
77 { }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82 {
83  return ( &particle == G4Positron::Positron() );
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 // Build cross section and mean free path tables
90 //here no tables, just calling PrintInfoDefinition
91 {
92  CurrentSigma = 0.0;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97 
99 // Set the factor to artificially increase the cross section
100 {
102  G4cout << "The cross section for AnnihiToMuPair is artificially "
103  << "increased by the CrossSecFactor=" << CrossSecFactor << G4endl;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
109 // Calculates the microscopic cross section in GEANT4 internal units.
110 // It gives a good description from threshold to 1000 GeV
111 {
112  static const G4double Mmuon = G4MuonPlus::MuonPlus()->GetPDGMass();
113  static const G4double Rmuon = elm_coupling/Mmuon; //classical particle radius
114  static const G4double Sig0 = pi*Rmuon*Rmuon/3.; //constant in crossSection
115 
116  G4double CrossSection = 0.;
117  if (Epos < LowestEnergyLimit) return CrossSection;
118 
119  G4double xi = LowestEnergyLimit/Epos;
120  G4double SigmaEl = Sig0*xi*(1.+xi/2.)*sqrt(1.-xi); // per electron
121  CrossSection = SigmaEl*Z; // number of electrons per atom
122  return CrossSection;
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
128  const G4Material* aMaterial)
129 {
130  const G4ElementVector* theElementVector = aMaterial->GetElementVector();
131  const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
132 
133  G4double SIGMA = 0.0;
134 
135  for ( size_t i=0 ; i < aMaterial->GetNumberOfElements() ; ++i )
136  {
137  G4double AtomicZ = (*theElementVector)[i]->GetZ();
138  SIGMA += NbOfAtomsPerVolume[i] *
139  ComputeCrossSectionPerAtom(PositronEnergy,AtomicZ);
140  }
141  return SIGMA;
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
148 
149 // returns the positron mean free path in GEANT4 internal units
150 
151 {
152  const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
153  G4double PositronEnergy = aDynamicPositron->GetKineticEnergy()
154  +electron_mass_c2;
155  G4Material* aMaterial = aTrack.GetMaterial();
156  CurrentSigma = CrossSectionPerVolume(PositronEnergy, aMaterial);
157 
158  // increase the CrossSection by CrossSecFactor (default 1)
159  G4double mfp = DBL_MAX;
161 
162  return mfp;
163 }
164 
165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166 
168  const G4Step& aStep)
169 //
170 // generation of e+e- -> mu+mu-
171 //
172 {
173 
174  aParticleChange.Initialize(aTrack);
175  static const G4double Mele=electron_mass_c2;
176  static const G4double Mmuon=G4MuonPlus::MuonPlus()->GetPDGMass();
177 
178  // current Positron energy and direction, return if energy too low
179  const G4DynamicParticle *aDynamicPositron = aTrack.GetDynamicParticle();
180  G4double Epos = aDynamicPositron->GetKineticEnergy() + Mele;
181 
182  // test of cross section
184  CrossSectionPerVolume(Epos, aTrack.GetMaterial()))
185  {
186  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
187  }
188 
189  if (Epos < LowestEnergyLimit) {
190  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
191  }
192 
193  G4ParticleMomentum PositronDirection =
194  aDynamicPositron->GetMomentumDirection();
195  G4double xi = LowestEnergyLimit/Epos; // xi is always less than 1,
196  // goes to 0 at high Epos
197 
198  // generate cost
199  //
200  G4double cost;
201  do { cost = 2.*G4UniformRand()-1.; }
202  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
203  while (2.*G4UniformRand() > 1.+xi+cost*cost*(1.-xi) );
204  //1+cost**2 at high Epos
205  G4double sint = sqrt(1.-cost*cost);
206 
207  // generate phi
208  //
209  G4double phi=2.*pi*G4UniformRand();
210 
211  G4double Ecm = sqrt(0.5*Mele*(Epos+Mele));
212  G4double Pcm = sqrt(Ecm*Ecm-Mmuon*Mmuon);
213  G4double beta = sqrt((Epos-Mele)/(Epos+Mele));
214  G4double gamma = Ecm/Mele; // =sqrt((Epos+Mele)/(2.*Mele));
215  G4double Pt = Pcm*sint;
216 
217  // energy and momentum of the muons in the Lab
218  //
219  G4double EmuPlus = gamma*( Ecm+cost*beta*Pcm);
220  G4double EmuMinus = gamma*( Ecm-cost*beta*Pcm);
221  G4double PmuPlusZ = gamma*(beta*Ecm+cost* Pcm);
222  G4double PmuMinusZ = gamma*(beta*Ecm-cost* Pcm);
223  G4double PmuPlusX = Pt*cos(phi);
224  G4double PmuPlusY = Pt*sin(phi);
225  G4double PmuMinusX =-Pt*cos(phi);
226  G4double PmuMinusY =-Pt*sin(phi);
227  // absolute momenta
228  G4double PmuPlus = sqrt(Pt*Pt+PmuPlusZ *PmuPlusZ );
229  G4double PmuMinus = sqrt(Pt*Pt+PmuMinusZ*PmuMinusZ);
230 
231  // mu+ mu- directions for Positron in z-direction
232  //
234  MuPlusDirection ( PmuPlusX/PmuPlus, PmuPlusY/PmuPlus, PmuPlusZ/PmuPlus );
236  MuMinusDirection(PmuMinusX/PmuMinus,PmuMinusY/PmuMinus,PmuMinusZ/PmuMinus);
237 
238  // rotate to actual Positron direction
239  //
240  MuPlusDirection.rotateUz(PositronDirection);
241  MuMinusDirection.rotateUz(PositronDirection);
242 
244  // create G4DynamicParticle object for the particle1
245  G4DynamicParticle* aParticle1= new G4DynamicParticle(
246  G4MuonPlus::MuonPlus(),MuPlusDirection,EmuPlus-Mmuon);
247  aParticleChange.AddSecondary(aParticle1);
248  // create G4DynamicParticle object for the particle2
249  G4DynamicParticle* aParticle2= new G4DynamicParticle(
250  G4MuonMinus::MuonMinus(),MuMinusDirection,EmuMinus-Mmuon);
251  aParticleChange.AddSecondary(aParticle2);
252 
253  // Kill the incident positron
254  //
257 
258  return &aParticleChange;
259 }
260 
261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
262 
264 {
265  G4String comments ="e+e->mu+mu- annihilation, atomic e- at rest, SubType=.";
266  G4cout << G4endl << GetProcessName() << ": " << comments
267  << GetProcessSubType() << G4endl;
268  G4cout << " threshold at " << LowestEnergyLimit/GeV << " GeV"
269  << " good description up to "
270  << HighestEnergyLimit/TeV << " TeV for all Z." << G4endl;
271 }
272 
273 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
const G4DynamicParticle * GetDynamicParticle() const
G4double ComputeCrossSectionPerAtom(G4double PositronEnergy, G4double AtomicZ)
G4bool IsApplicable(const G4ParticleDefinition &)
G4AnnihiToMuPair(const G4String &processName="AnnihiToMuPair", G4ProcessType type=fElectromagnetic)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double GeV
Definition: G4SIunits.hh:214
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
virtual void Initialize(const G4Track &)
static G4Positron * Positron()
Definition: G4Positron.cc:94
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *)
G4double GetPDGMass() const
static const double pi
Definition: G4SIunits.hh:74
void SetNumberOfSecondaries(G4int totSecondaries)
void ProposeEnergy(G4double finalEnergy)
#define DBL_MIN
Definition: templates.hh:75
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static const G4double fac
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4double CrossSectionPerVolume(G4double PositronEnergy, const G4Material *)
void SetCrossSecFactor(G4double fac)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4ThreeVector G4ParticleMomentum
#define DBL_MAX
Definition: templates.hh:83
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ProcessType