Geant4  10.02.p01
G4MuonDecayChannelWithSpin.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 // GEANT 4 class header file
28 //
29 // History:
30 // 17 August 2004 P.Gumplinger and T.MacPhail
31 // samples Michel spectrum including 1st order
32 // radiative corrections
33 // Reference: Florian Scheck "Muon Physics", in Physics Reports
34 // (Review Section of Physics Letters) 44, No. 4 (1978)
35 // 187-248. North-Holland Publishing Company, Amsterdam
36 // at page 210 cc.
37 //
38 // W.E. Fisher and F. Scheck, Nucl. Phys. B83 (1974) 25.
39 //
40 // ------------------------------------------------------------
41 //
43 
44 #include "G4PhysicalConstants.hh"
45 #include "G4SystemOfUnits.hh"
46 #include "Randomize.hh"
47 
48 #include "G4DecayProducts.hh"
49 #include "G4LorentzVector.hh"
50 
53  EMMU( 0.*MeV),
54  EMASS( 0.*MeV)
55 {
56 }
57 
59  G4double theBR)
60  : G4MuonDecayChannel(theParentName,theBR),
61  EMMU( 0.*MeV),
62  EMASS( 0.*MeV)
63 {
64 }
65 
67 {
68 }
69 
71  G4MuonDecayChannel(right)
72 {
73  EMMU = right.EMMU;
74  EMASS = right.EMASS;
75 }
76 
78 {
79  if (this != &right) {
81  verboseLevel = right.verboseLevel;
82  rbranch = right.rbranch;
83 
84  // copy parent name
85  parent_name = new G4String(*right.parent_name);
86 
87  // clear daughters_name array
89 
90  // recreate array
92  if ( numberOfDaughters >0 ) {
95  //copy daughters name
96  for (G4int index=0; index < numberOfDaughters; index++) {
97  daughters_name[index] = new G4String(*right.daughters_name[index]);
98  }
99  }
100  EMMU = right.EMMU;
101  EMASS = right.EMASS;
102  }
103  return *this;
104 }
105 
106 
108 {
109  // This version assumes V-A coupling with 1st order radiative correctons,
110  // the standard model Michel parameter values, but
111  // gives incorrect energy spectrum for neutrinos
112 
113 #ifdef G4VERBOSE
114  if (GetVerboseLevel()>1) G4cout << "G4MuonDecayChannelWithSpin::DecayIt ";
115 #endif
116 
117  if (G4MT_parent == 0) FillParent();
118  if (G4MT_daughters == 0) FillDaughters();
119 
120  // parent mass
121  G4double parentmass = G4MT_parent->GetPDGMass();
122 
123  EMMU = parentmass;
124 
125  //daughters'mass
126  G4double daughtermass[3];
127  G4double sumofdaughtermass = 0.0;
128  for (G4int index=0; index<3; index++){
129  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
130  sumofdaughtermass += daughtermass[index];
131  }
132 
133  EMASS = daughtermass[0];
134 
135  //create parent G4DynamicParticle at rest
136  G4ThreeVector dummy;
137  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
138  //create G4Decayproducts
139  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
140  delete parentparticle;
141 
142  // calcurate electron energy
143 
144  G4double michel_rho = 0.75; //Standard Model Michel rho
145  G4double michel_delta = 0.75; //Standard Model Michel delta
146  G4double michel_xsi = 1.00; //Standard Model Michel xsi
147  G4double michel_eta = 0.00; //Standard Model eta
148 
149  G4double rndm, x, ctheta;
150 
151  G4double FG;
152  G4double FG_max = 2.00;
153 
154  G4double W_mue = (EMMU*EMMU+EMASS*EMASS)/(2.*EMMU);
155  G4double x0 = EMASS/W_mue;
156 
157  G4double x0_squared = x0*x0;
158 
159  // ***************************************************
160  // x0 <= x <= 1. and -1 <= y <= 1
161  //
162  // F(x,y) = f(x)*g(x,y); g(x,y) = 1.+g(x)*y
163  // ***************************************************
164 
165  // ***** sampling F(x,y) directly (brute force) *****
166 
167  const size_t MAX_LOOP=10000;
168  for (size_t loop_count =0; loop_count <MAX_LOOP; ++loop_count){
169 
170  // Sample the positron energy by sampling from F
171 
172  rndm = G4UniformRand();
173 
174  x = x0 + rndm*(1.-x0);
175 
176  G4double x_squared = x*x;
177 
178  G4double F_IS, F_AS, G_IS, G_AS;
179 
180  F_IS = 1./6.*(-2.*x_squared+3.*x-x0_squared);
181  F_AS = 1./6.*std::sqrt(x_squared-x0_squared)*(2.*x-2.+std::sqrt(1.-x0_squared));
182 
183  G_IS = 2./9.*(michel_rho-0.75)*(4.*x_squared-3.*x-x0_squared);
184  G_IS = G_IS + michel_eta*(1.-x)*x0;
185 
186  G_AS = 3.*(michel_xsi-1.)*(1.-x);
187  G_AS = G_AS+2.*(michel_xsi*michel_delta-0.75)*(4.*x-4.+std::sqrt(1.-x0_squared));
188  G_AS = 1./9.*std::sqrt(x_squared-x0_squared)*G_AS;
189 
190  F_IS = F_IS + G_IS;
191  F_AS = F_AS + G_AS;
192 
193  // *** Radiative Corrections ***
194 
195  G4double R_IS = F_c(x,x0);
196 
197  G4double F = 6.*F_IS + R_IS/std::sqrt(x_squared-x0_squared);
198 
199  // *** Radiative Corrections ***
200 
201  G4double R_AS = F_theta(x,x0);
202 
203  rndm = G4UniformRand();
204 
205  ctheta = 2.*rndm-1.;
206 
207  G4double G = 6.*F_AS - R_AS/std::sqrt(x_squared-x0_squared);
208 
209  FG = std::sqrt(x_squared-x0_squared)*F*(1.+(G/F)*ctheta);
210 
211  if(FG>FG_max){
212  G4cout<<"***Problem in Muon Decay *** : FG > FG_max"<<G4endl;
213  FG_max = FG;
214  }
215 
216  rndm = G4UniformRand();
217 
218  if (FG >= rndm*FG_max) break;
219  }
220 
221  G4double energy = x * W_mue;
222 
223  rndm = G4UniformRand();
224 
225  G4double phi = twopi * rndm;
226 
227  if(energy < EMASS) energy = EMASS;
228 
229  // calculate daughter momentum
230  G4double daughtermomentum[3];
231 
232  daughtermomentum[0] = std::sqrt(energy*energy - EMASS*EMASS);
233 
234  G4double stheta = std::sqrt(1.-ctheta*ctheta);
235  G4double cphi = std::cos(phi);
236  G4double sphi = std::sin(phi);
237 
238  //Coordinates of the decay positron with respect to the muon spin
239 
240  G4double px = stheta*cphi;
241  G4double py = stheta*sphi;
242  G4double pz = ctheta;
243 
244  G4ThreeVector direction0(px,py,pz);
245 
246  direction0.rotateUz(parent_polarization);
247 
248  G4DynamicParticle * daughterparticle0
249  = new G4DynamicParticle( G4MT_daughters[0], daughtermomentum[0]*direction0);
250 
251  products->PushProducts(daughterparticle0);
252 
253 
254  // daughter 1 ,2 (neutrinos)
255  // create neutrinos in the C.M frame of two neutrinos
256  G4double energy2 = parentmass*(1.0 - x/2.0);
257  G4double vmass = std::sqrt((energy2-daughtermomentum[0])*(energy2+daughtermomentum[0]));
258  G4double beta = -1.0*daughtermomentum[0]/energy2;
259  G4double costhetan = 2.*G4UniformRand()-1.0;
260  G4double sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
261  G4double phin = twopi*G4UniformRand()*rad;
262  G4double sinphin = std::sin(phin);
263  G4double cosphin = std::cos(phin);
264 
265  G4ThreeVector direction1(sinthetan*cosphin,sinthetan*sinphin,costhetan);
266  G4DynamicParticle * daughterparticle1
267  = new G4DynamicParticle( G4MT_daughters[1], direction1*(vmass/2.));
268  G4DynamicParticle * daughterparticle2
269  = new G4DynamicParticle( G4MT_daughters[2], direction1*(-1.0*vmass/2.));
270 
271  // boost to the muon rest frame
272  G4LorentzVector p4;
273  p4 = daughterparticle1->Get4Momentum();
274  p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
275  daughterparticle1->Set4Momentum(p4);
276  p4 = daughterparticle2->Get4Momentum();
277  p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
278  daughterparticle2->Set4Momentum(p4);
279  products->PushProducts(daughterparticle1);
280  products->PushProducts(daughterparticle2);
281  daughtermomentum[1] = daughterparticle1->GetTotalMomentum();
282  daughtermomentum[2] = daughterparticle2->GetTotalMomentum();
283 
284  // output message
285 #ifdef G4VERBOSE
286  if (GetVerboseLevel()>1) {
287  G4cout << "G4MuonDecayChannelWithSpin::DecayIt ";
288  G4cout << " create decay products in rest frame " <<G4endl;
289  products->DumpInfo();
290  }
291 #endif
292  return products;
293 }
294 
296 
297  G4int n_max = (int)(100.*x);
298 
299  if(n_max<10)n_max=10;
300 
301  G4double L2 = 0.0;
302 
303  for(G4int n=1; n<=n_max; n++){
304  L2 += std::pow(x,n)/(n*n);
305  }
306 
307  G4double omega = std::log(EMMU/EMASS);
308 
309  G4double r_c;
310 
311  r_c = 2.*L2-(pi*pi/3.)-2.;
312  r_c = r_c + omega * (1.5+2.*std::log((1.-x)/x));
313  r_c = r_c - std::log(x)*(2.*std::log(x)-1.);
314  r_c = r_c + (3.*std::log(x)-1.-1./x)*std::log(1.-x);
315 
316  return r_c;
317 }
static const double MeV
Definition: G4SIunits.hh:211
CLHEP::Hep3Vector G4ThreeVector
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition ** G4MT_daughters
G4double F_theta(G4double x, G4double x0)
int G4int
Definition: G4Types.hh:78
virtual G4DecayProducts * DecayIt(G4double)
G4MuonDecayChannelWithSpin & operator=(const G4MuonDecayChannelWithSpin &)
G4double GetTotalMomentum() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4String kinematics_name
void DumpInfo() const
static const double twopi
Definition: G4SIunits.hh:75
G4String * parent_name
const G4int n
G4ThreeVector parent_polarization
G4LorentzVector Get4Momentum() const
static const double rad
Definition: G4SIunits.hh:148
void Set4Momentum(const G4LorentzVector &momentum)
G4int GetVerboseLevel() const
G4double GetPDGMass() const
static const double pi
Definition: G4SIunits.hh:74
G4double F_c(G4double x, G4double x0)
G4double energy(const ThreeVector &p, const G4double m)
const G4double x[NPOINTSGL]
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4String ** daughters_name
CLHEP::HepLorentzVector G4LorentzVector