Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4MuMinusCaptureCascade Class Reference

#include <G4MuMinusCaptureCascade.hh>

Public Member Functions

 G4MuMinusCaptureCascade ()
 
 ~G4MuMinusCaptureCascade ()
 
G4int DoCascade (const G4double Z, const G4double A, G4GHEKinematicsVector *Cascade)
 
void DoBoundMuonMinusDecay (G4double Z, G4int *nCascade, G4GHEKinematicsVector *Cascade)
 
G4double GetKShellEnergy (G4double Z)
 
G4ThreeVectorGetRandomVec ()
 

Detailed Description

Definition at line 48 of file G4MuMinusCaptureCascade.hh.

Constructor & Destructor Documentation

G4MuMinusCaptureCascade::G4MuMinusCaptureCascade ( )

Definition at line 55 of file G4MuMinusCaptureCascade.cc.

References G4Electron::Electron(), G4Gamma::Gamma(), G4ParticleDefinition::GetPDGMass(), and G4MuonMinus::MuonMinus().

56 {
57  theElectron = G4Electron::Electron();
58  theGamma = G4Gamma::Gamma();
59  Emass = theElectron->GetPDGMass();
60  MuMass = G4MuonMinus::MuonMinus()->GetPDGMass();
61 }
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetPDGMass() const
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4MuMinusCaptureCascade::~G4MuMinusCaptureCascade ( )

Definition at line 65 of file G4MuMinusCaptureCascade.cc.

66 { }

Member Function Documentation

void G4MuMinusCaptureCascade::DoBoundMuonMinusDecay ( G4double  Z,
G4int nCascade,
G4GHEKinematicsVector Cascade 
)

Definition at line 181 of file G4MuMinusCaptureCascade.cc.

References G4AntiNeutrinoE::AntiNeutrinoE(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), CLHEP::HepLorentzVector::e(), G4UniformRand, GetKShellEnergy(), GetRandomVec(), CLHEP::HepLorentzVector::mag2(), G4NeutrinoMu::NeutrinoMu(), CLHEP::Hep3Vector::unit(), CLHEP::HepLorentzVector::vect(), and test::x.

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

184 {
185  // Simulation on Decay of mu- on a K-shell of the muonic atom
186  G4double xmax = ( 1.0 + Emass*Emass/ (MuMass*MuMass) );
187  G4double xmin = 2.0*Emass/MuMass;
188  G4double KEnergy = GetKShellEnergy(Z);
189  /*
190  G4cout << "G4MuMinusCaptureCascade::DoBoundMuonMinusDecay"
191  << " XMAX= " << xmax
192  << " Ebound= " << KEnergy
193  << G4endl;
194  */
195  G4double pmu = std::sqrt(KEnergy*(KEnergy + 2.0*MuMass));
196  G4double emu = KEnergy + MuMass;
197  G4ThreeVector moment = GetRandomVec();
198  G4LorentzVector MU(pmu*moment,emu);
199  G4ThreeVector bst = MU.boostVector();
200 
201  G4double Eelect, Pelect, x, ecm;
202  G4LorentzVector EL, NN;
203  // Calculate electron energy
204  do {
205  do {
206  x = xmin + (xmax-xmin)*G4UniformRand();
207  } while (G4UniformRand() > (3.0 - 2.0*x)*x*x );
208  Eelect = x*MuMass*0.5;
209  Pelect = 0.0;
210  if(Eelect > Emass) {
211  Pelect = std::sqrt( Eelect*Eelect - Emass*Emass );
212  } else {
213  Pelect = 0.0;
214  Eelect = Emass;
215  }
216  G4ThreeVector e_mom = GetRandomVec();
217  EL = G4LorentzVector(Pelect*e_mom,Eelect);
218  EL.boost(bst);
219  Eelect = EL.e() - Emass - 2.0*KEnergy;
220  //
221  // Calculate rest frame parameters of 2 neutrinos
222  //
223  NN = MU - EL;
224  ecm = NN.mag2();
225  } while (Eelect < 0.0 || ecm < 0.0);
226 
227  //
228  // Create electron
229  //
230  moment = std::sqrt(Eelect * (Eelect + 2.0*Emass))*(EL.vect().unit());
231  AddNewParticle(theElectron, moment, Emass, nCascade, Cascade);
232  //
233  // Create Neutrinos
234  //
235  ecm = 0.5*std::sqrt(ecm);
236  bst = NN.boostVector();
237  G4ThreeVector p1 = ecm * GetRandomVec();
238  G4LorentzVector N1 = G4LorentzVector(p1,ecm);
239  N1.boost(bst);
240  G4ThreeVector p1lab = N1.vect();
241  AddNewParticle(G4AntiNeutrinoE::AntiNeutrinoE(),p1lab,0.0,nCascade,Cascade);
242  NN -= N1;
243  G4ThreeVector p2lab = NN.vect();
244  AddNewParticle(G4NeutrinoMu::NeutrinoMu(),p2lab,0.0,nCascade,Cascade);
245 
246  return;
247 }
Hep3Vector boostVector() const
static G4AntiNeutrinoE * AntiNeutrinoE()
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
static G4NeutrinoMu * NeutrinoMu()
Definition: G4NeutrinoMu.cc:85
HepLorentzVector & boost(double, double, double)
G4double GetKShellEnergy(G4double Z)
double mag2() const
double G4double
Definition: G4Types.hh:76
CLHEP::HepLorentzVector G4LorentzVector
G4int G4MuMinusCaptureCascade::DoCascade ( const G4double  Z,
const G4double  A,
G4GHEKinematicsVector Cascade 
)

Definition at line 114 of file G4MuMinusCaptureCascade.cc.

References python.hepunit::electron_mass_c2, python.hepunit::eV, G4UniformRand, GetKShellEnergy(), GetRandomVec(), and var().

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

116 {
117  // Inicialization - cascade start from 14th level
118  // N.C.Mukhopadhyay Phy. Rep. 30 (1977) 1.
119  G4int nPart = 0;
120  G4double EnergyLevel[14];
121 
122  G4double mass = MuMass * massA / (MuMass + massA) ;
123 
124  const G4double KEnergy = 13.6 * eV * Z * Z * mass/ electron_mass_c2;
125 
126  EnergyLevel[0] = GetKShellEnergy(Z);
127  for( G4int i = 2; i < 15; i++ ) {
128  EnergyLevel[i-1] = KEnergy / (i*i) ;
129  }
130 
131  G4int nElec = G4int(Z);
132  G4int nAuger = 1;
133  G4int nLevel = 13;
134  G4double DeltaE;
135  G4double pGamma = Z*Z*Z*Z;
136 
137  // Capture on 14-th level
138  G4double ptot = std::sqrt(EnergyLevel[13]*(EnergyLevel[13] + 2.0*Emass));
139  G4ThreeVector moment = ptot * GetRandomVec();
140 
141  AddNewParticle(theElectron,moment,Emass,&nPart,Cascade);
142 
143  // Emit new photon or electron
144  // Simplified model for probabilities
145  // N.C.Mukhopadhyay Phy. Rep. 30 (1977) 1.
146  do {
147 
148  // case of Auger electrons
149  if((nAuger < nElec) && ((pGamma + 10000.0) * G4UniformRand() < 10000.0) ) {
150  nAuger++;
151  DeltaE = EnergyLevel[nLevel-1] - EnergyLevel[nLevel];
152  nLevel--;
153 
154  ptot = std::sqrt(DeltaE * (DeltaE + 2.0*Emass));
155  moment = ptot * GetRandomVec();
156 
157  AddNewParticle(theElectron, moment, Emass, &nPart, Cascade);
158 
159  } else {
160 
161  // Case of photon cascade, probabilities from
162  // C.S.Wu and L.Wilets, Ann. Rev. Nuclear Sci. 19 (1969) 527.
163 
164  G4double var = (10.0 + G4double(nLevel - 1) ) * G4UniformRand();
165  G4int iLevel = nLevel - 1 ;
166  if(var > 10.0) iLevel -= G4int(var-10.0) + 1;
167  if( iLevel < 0 ) iLevel = 0;
168  DeltaE = EnergyLevel[iLevel] - EnergyLevel[nLevel];
169  nLevel = iLevel;
170  moment = DeltaE * GetRandomVec();
171  AddNewParticle(theGamma, moment, 0.0, &nPart, Cascade);
172  }
173 
174  } while( nLevel > 0 );
175 
176  return nPart;
177 }
real *8 function var(A, B, C, D)
Definition: dpm25nuc1.f:4649
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
float electron_mass_c2
Definition: hepunit.py:274
G4double GetKShellEnergy(G4double Z)
double G4double
Definition: G4Types.hh:76
G4double G4MuMinusCaptureCascade::GetKShellEnergy ( G4double  Z)

Definition at line 70 of file G4MuMinusCaptureCascade.cc.

Referenced by DoBoundMuonMinusDecay(), and DoCascade().

71 {
72  // Calculate the Energy of K Mesoatom Level for this Element using
73  // the Energy of Hydrogen Atom taken into account finite size of the
74  // nucleus (V.Ivanchenko)
75  const G4int ListK = 28;
76  const G4double ListZK[ListK] = {
77  1., 2., 4., 6., 8., 11., 14., 17., 18., 21., 24.,
78  26., 29., 32., 38., 40., 41., 44., 49., 53., 55.,
79  60., 65., 70., 75., 81., 85., 92.};
80  const G4double ListKEnergy[ListK] = {
81  0.00275, 0.011, 0.043, 0.098, 0.173, 0.326,
82  0.524, 0.765, 0.853, 1.146, 1.472,
83  1.708, 2.081, 2.475, 3.323, 3.627,
84  3.779, 4.237, 5.016, 5.647, 5.966,
85  6.793, 7.602, 8.421, 9.249, 10.222,
86  10.923,11.984};
87 
88  // Energy with finit size corrections
89  G4double KEnergy = GetLinApprox(ListK,ListZK,ListKEnergy,Z);
90 
91  return KEnergy;
92 }
int G4int
Definition: G4Types.hh:78
double G4double
Definition: G4Types.hh:76
G4ThreeVector & G4MuMinusCaptureCascade::GetRandomVec ( )
inline

Definition at line 111 of file G4MuMinusCaptureCascade.hh.

References G4UniformRand.

Referenced by DoBoundMuonMinusDecay(), and DoCascade().

112 {
113  //
114  // generate uniform vector
115  //
116  G4double cost = 2.0 * G4UniformRand() - 1.0;
117  G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
118  G4double Phi = CLHEP::twopi * G4UniformRand();
119  randomVect = G4ThreeVector(sint * std::cos(Phi), sint * std::sin(Phi), cost);
120  return randomVect;
121 }
CLHEP::Hep3Vector G4ThreeVector
#define G4UniformRand()
Definition: Randomize.hh:87
double G4double
Definition: G4Types.hh:76

The documentation for this class was generated from the following files: