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