Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuonMinusBoundDecay.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// $Id$
27//
28//-----------------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32// File name: G4MuonMinusBoundDecay
33//
34// Author: V.Ivanchenko ([email protected])
35//
36// Creation date: 24 April 2012 on base of G4MuMinusCaptureAtRest
37//
38// Modified:
39//
40//----------------------------------------------------------------------
41
43#include "Randomize.hh"
44#include "G4RandomDirection.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4ThreeVector.hh"
48#include "G4MuonMinus.hh"
49#include "G4Electron.hh"
50#include "G4NeutrinoMu.hh"
51#include "G4AntiNeutrinoE.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
56 : G4HadronicInteraction("muMinusBoundDeacy")
57{
58 fMuMass = G4MuonMinus::MuonMinus()->GetPDGMass();
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
64{}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
70 G4Nucleus& targetNucleus)
71{
72 result.Clear();
73 G4int Z = targetNucleus.GetZ_asInt();
74 G4int A = targetNucleus.GetA_asInt();
75
76 // Decide on Decay or Capture, and doit.
77 G4double lambdac = GetMuonCaptureRate(Z, A);
78 G4double lambdad = GetMuonDecayRate(Z);
79 G4double lambda = lambdac + lambdad;
80
81 // === sample capture time and change time of projectile
82
83 G4double time = -std::log(G4UniformRand()) / lambda;
84 G4HadProjectile* p = const_cast<G4HadProjectile*>(&projectile);
85 p->SetGlobalTime(time);
86
87 //G4cout << "lambda= " << lambda << " lambdac= " << lambdac
88 //<< " t= " << time << G4endl;
89
90 // cascade
91 if( G4UniformRand()*lambda < lambdac) {
93
94 } else {
95
96 // Simulation on Decay of mu- on a K-shell of the muonic atom
98 G4double xmax = 1 + electron_mass_c2*electron_mass_c2/(fMuMass*fMuMass);
99 G4double xmin = 2.0*electron_mass_c2/fMuMass;
100 G4double KEnergy = projectile.GetBoundEnergy();
101
102 /*
103 G4cout << "G4MuonMinusBoundDecay::ApplyYourself"
104 << " XMAX= " << xmax << " Ebound= " << KEnergy<< G4endl;
105 */
106 G4double pmu = std::sqrt(KEnergy*(KEnergy + 2.0*fMuMass));
107 G4double emu = KEnergy + fMuMass;
109 G4LorentzVector MU(pmu*dir, emu);
110 G4ThreeVector bst = MU.boostVector();
111
112 G4double Eelect, Pelect, x, ecm;
113 G4LorentzVector EL, NN;
114 // Calculate electron energy
115 do {
116 do {
117 x = xmin + (xmax-xmin)*G4UniformRand();
118 } while (G4UniformRand() > (3.0 - 2.0*x)*x*x );
119 Eelect = x*fMuMass*0.5;
120 Pelect = 0.0;
121 if(Eelect > electron_mass_c2) {
122 Pelect = std::sqrt(Eelect*Eelect - electron_mass_c2*electron_mass_c2);
123 } else {
124 Pelect = 0.0;
125 Eelect = electron_mass_c2;
126 }
127 dir = G4RandomDirection();
128 EL = G4LorentzVector(Pelect*dir,Eelect);
129 EL.boost(bst);
130 Eelect = EL.e() - electron_mass_c2 - 2.0*KEnergy;
131 //
132 // Calculate rest frame parameters of 2 neutrinos
133 //
134 NN = MU - EL;
135 ecm = NN.mag2();
136 } while (Eelect < 0.0 || ecm < 0.0);
137
138 //
139 // Create electron
140 //
142 EL.vect().unit(),
143 Eelect);
144
145 AddNewParticle(dp, time);
146 //
147 // Create Neutrinos
148 //
149 ecm = 0.5*std::sqrt(ecm);
150 bst = NN.boostVector();
151 G4ThreeVector p1 = ecm * G4RandomDirection();
152 G4LorentzVector N1 = G4LorentzVector(p1,ecm);
153 N1.boost(bst);
155 AddNewParticle(dp, time);
156 NN -= N1;
158 AddNewParticle(dp, time);
159 }
160 return &result;
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
165G4double G4MuonMinusBoundDecay::GetMuonCaptureRate(G4int Z, G4int A)
166{
167 // Initialized data
168 static G4double zeff[101] = { 0.,
169 1.,1.98,2.95,3.89,4.8,5.72,6.61,7.49,8.32,9.12,9.95,10.69,11.48,12.22,
170 12.91,13.64,14.24,14.89,15.53,16.15,16.75,17.38,18.04,18.49,
171 19.06,19.59,20.1,20.66,21.12,21.61,22.02,22.43,22.84,23.24,
172 23.65,24.06,24.47,24.85,25.23,25.61,25.99,26.37,26.69,27.,
173 27.32,27.63,27.95,28.2,28.42,28.64,28.79,29.03,29.27,29.51,
174 29.75,29.99,30.2,30.36,30.53,30.69,30.85,31.01,31.18,31.34,
175 31.48,31.62,31.76,31.9,32.05,32.19,32.33,32.47,32.61,32.76,
176 32.94,33.11,33.29,33.46,33.64,33.81,34.21,34.18,34.,34.1,
177 34.21,34.31,34.42,34.52,34.63,34.73,34.84,34.94,35.04,35.15,
178 35.25,35.36,35.46,35.57,35.67,35.78 };
179
180 // Mu- capture data from B.B.Balashov, G.Ya.Korenman, P.A.Eramgan
181 // Atomizdat, 1978. (Experimental capture velocities)
182 // Data for Hydrogen from Phys. Rev. Lett. 99(2007)032002
183 // Data for Helium from Phys. Rep. 354(2001)243
184
185 const size_t ListZE = 67;
186 static G4int ListZExp[ListZE] = { 1, 2,
187 3, 4, 5, 6, 7, 8, 9, 10, 11, 12,
188 13, 14, 15, 16, 17, 18, 19, 20, 22, 23,
189 24, 25, 26, 27, 28, 31, 32, 33, 34, 37,
190 38, 39, 40, 41, 42, 45, 46, 47, 48, 49,
191 50, 51, 52, 53, 55, 56, 57, 58, 59, 60,
192 62, 64, 65, 67, 72, 73, 74, 80, 81, 82,
193 83, 90, 92, 93};
194
195 static G4double ListCaptureVel[ListZE] = { 0.000725, 0.000356,
196 0.0057, 0.010, 0.0258, 0.0371, 0.0644,
197 0.0974, 0.144, 0.250, 0.386, 0.479,
198 0.700, 0.849, 1.119, 1.338, 1.40,
199 1.30, 1.98, 2.45, 2.60, 3.19,
200 3.29, 3.91, 4.41, 4.96, 5.74,
201 5.68, 5.53, 6.06, 5.69, 6.89,
202 7.25, 7.89, 8.59, 10.40, 9.22,
203 10.01, 10.00, 10.88, 10.62, 11.37,
204 10.68, 10.49, 9.06, 11.20, 10.98,
205 10.18, 10.71, 11.44, 13.45, 12.32,
206 12.22, 12.09, 12.73, 12.95, 13.03,
207 12.86, 13.13, 13.39, 12.74, 13.78,
208 13.02, 13.26, 13.10, 14.00, 14.70};
209
210
211 // Local variables
212 G4double zeff2, xmu, a2ze, r1, r2;
214
215 // == Effective charges from Ford and Wills Nucl Phys 35(1962)295.
216 // == Untabulated charges are interpolated.
217 // == Mu capture lifetime (Goulard and Primakoff PRC10(1974)2034.
218
219 G4int i = Z;
220 if(i > 100) { i = 100; }
221
222 const G4double b0a = -.03;
223 const G4double b0b = -.25;
224 const G4double b0c = 3.24;
225 const G4double t1 = 875.e-10;
226 r1 = zeff[i];
227 zeff2 = r1 * r1;
228
229 // ^-4 -> ^-5 suggested by user
230 xmu = zeff2 * 2.663e-5;
231 a2ze = 0.5 * A / Z;
232 r2 = 1.0 - xmu;
233 lambda = t1 * zeff2 * zeff2 * (r2 * r2) * (1.0 - (1.0 - xmu) * .75704) *
234 (a2ze * b0a + 1.0 - (a2ze - 1.0) * b0b -
235 (2 * (A - Z) + std::fabs(a2ze - 1.) ) * b0c / G4double(A * 4) );
236
237 // == Mu capture data are taken if exist
238 for (size_t j = 0; j < ListZE; ++j) {
239 if( ListZExp[j] == i + 1) {
240 lambda = ListCaptureVel[j] / microsecond;
241 break;
242 }
243 }
244
245 return lambda;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249
250G4double G4MuonMinusBoundDecay::GetMuonDecayRate(G4int Z)
251{
252 // Decay time on K-shell
253 // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
254
255 G4double lambda = 1.0;
256 if(Z > 1) {
257 G4double x = Z*fine_structure_const;
258 lambda -= 2.5 * x * x;
259 if( 0.5 > lambda ) { lambda = 0.5; }
260 } else {
261
262 // Published value 0.455851 - Phys. Rev. Lett. 99(2007)032002
263 lambda = 1.00151;
264 }
265 return lambda * 0.445164 / microsecond;
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
269
270void G4MuonMinusBoundDecay::ModelDescription(std::ostream& outFile) const
271{
272 outFile << "Sample probabilities of mu- nuclear capture of decay"
273 << " from K-shell orbit.\n"
274 << " Time of projectile is changed taking into account life time"
275 << " of muonic atom.\n"
276 << " If decay is sampled primary state become stopAndKill,"
277 << " else - isAlive.\n"
278 << "Based of reviews:\n"
279 << " N.C.Mukhopadhyay Phy. Rep. 30 (1977) 1.\n"
280 << " B.B.Balashov, G.Ya.Korenman, P.A.Eramgan, Atomizdat, 1978.\n";
281
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
285
@ isAlive
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
static G4AntiNeutrinoE * AntiNeutrinoE()
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetStatusChange(G4HadFinalStateStatus aS)
G4double GetBoundEnergy() const
void SetGlobalTime(G4double t)
void ModelDescription(std::ostream &outFile) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4NeutrinoMu * NeutrinoMu()
Definition: G4NeutrinoMu.cc:85
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115