source: trunk/source/particles/management/src/G4MuonDecayChannelWithSpin.cc

Last change on this file was 1340, checked in by garnier, 14 years ago

update ti head

File size: 8.5 KB
Line 
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//
42#include "G4MuonDecayChannelWithSpin.hh"
43
44#include "Randomize.hh"
45
46#include "G4DecayProducts.hh"
47#include "G4LorentzVector.hh"
48
49G4MuonDecayChannelWithSpin::G4MuonDecayChannelWithSpin(const G4String& theParentName, 
50                                                       G4double        theBR)
51                           : G4MuonDecayChannel(theParentName,theBR)
52{
53  EMMU  = 0.*MeV;
54  EMASS = 0.*MeV;
55}
56
57G4MuonDecayChannelWithSpin::~G4MuonDecayChannelWithSpin()
58{
59}
60
61G4DecayProducts *G4MuonDecayChannelWithSpin::DecayIt(G4double) 
62{
63  // This version assumes V-A coupling with 1st order radiative correctons,
64  //              the standard model Michel parameter values, but
65  //              gives incorrect energy spectrum for neutrinos
66
67#ifdef G4VERBOSE
68  if (GetVerboseLevel()>1) G4cout << "G4MuonDecayChannelWithSpin::DecayIt ";
69#endif
70
71  if (parent == 0) FillParent(); 
72  if (daughters == 0) FillDaughters();
73
74  // parent mass
75  G4double parentmass = parent->GetPDGMass();
76
77  EMMU = parentmass;
78
79  //daughters'mass
80  G4double daughtermass[3]; 
81  G4double sumofdaughtermass = 0.0;
82  for (G4int index=0; index<3; index++){
83    daughtermass[index] = daughters[index]->GetPDGMass();
84    sumofdaughtermass += daughtermass[index];
85  }
86
87  EMASS = daughtermass[0];
88
89  //create parent G4DynamicParticle at rest
90  G4ThreeVector dummy;
91  G4DynamicParticle * parentparticle = new G4DynamicParticle( parent, dummy, 0.0);
92  //create G4Decayproducts
93  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
94  delete parentparticle;
95
96  // calcurate electron energy
97
98  G4double michel_rho   = 0.75; //Standard Model Michel rho
99  G4double michel_delta = 0.75; //Standard Model Michel delta
100  G4double michel_xsi   = 1.00; //Standard Model Michel xsi
101  G4double michel_eta   = 0.00; //Standard Model eta
102
103  G4double rndm, x, ctheta;
104
105  G4double FG; 
106  G4double FG_max = 2.00;
107
108  G4double W_mue  = (EMMU*EMMU+EMASS*EMASS)/(2.*EMMU);
109  G4double x0     =           EMASS/W_mue;
110
111  G4double x0_squared = x0*x0;
112
113  // ***************************************************
114  //     x0 <= x <= 1.   and   -1 <= y <= 1
115  //
116  //     F(x,y) = f(x)*g(x,y);   g(x,y) = 1.+g(x)*y
117  // ***************************************************
118
119  // ***** sampling F(x,y) directly (brute force) *****
120
121  do{
122
123    // Sample the positron energy by sampling from F
124
125    rndm = G4UniformRand();
126
127    x = x0 + rndm*(1.-x0);
128
129    G4double x_squared = x*x;
130
131    G4double F_IS, F_AS, G_IS, G_AS;
132
133    F_IS = 1./6.*(-2.*x_squared+3.*x-x0_squared);
134    F_AS = 1./6.*std::sqrt(x_squared-x0_squared)*(2.*x-2.+std::sqrt(1.-x0_squared));
135
136    G_IS = 2./9.*(michel_rho-0.75)*(4.*x_squared-3.*x-x0_squared);
137    G_IS = G_IS + michel_eta*(1.-x)*x0;
138
139    G_AS = 3.*(michel_xsi-1.)*(1.-x);
140    G_AS = G_AS+2.*(michel_xsi*michel_delta-0.75)*(4.*x-4.+std::sqrt(1.-x0_squared));
141    G_AS = 1./9.*std::sqrt(x_squared-x0_squared)*G_AS;
142
143    F_IS = F_IS + G_IS;
144    F_AS = F_AS + G_AS;
145
146    // *** Radiative Corrections ***
147
148    G4double R_IS = F_c(x,x0);
149
150    G4double F = 6.*F_IS + R_IS/std::sqrt(x_squared-x0_squared);
151
152    // *** Radiative Corrections ***
153
154    G4double R_AS = F_theta(x,x0);
155
156    rndm = G4UniformRand();
157
158    ctheta = 2.*rndm-1.;
159
160    G4double G = 6.*F_AS - R_AS/std::sqrt(x_squared-x0_squared);
161
162    FG = std::sqrt(x_squared-x0_squared)*F*(1.+(G/F)*ctheta);
163
164    if(FG>FG_max){
165      G4cout<<"***Problem in Muon Decay *** : FG > FG_max"<<G4endl;
166      FG_max = FG;
167    }
168
169    rndm = G4UniformRand();
170
171  }while(FG<rndm*FG_max);
172
173  G4double energy = x * W_mue;
174
175  rndm = G4UniformRand();
176
177  G4double phi = twopi * rndm;
178
179  if(energy < EMASS) energy = EMASS;
180
181  // calculate daughter momentum
182  G4double daughtermomentum[3];
183
184  daughtermomentum[0] = std::sqrt(energy*energy - EMASS*EMASS);
185
186  G4double stheta = std::sqrt(1.-ctheta*ctheta);
187  G4double cphi = std::cos(phi);
188  G4double sphi = std::sin(phi);
189
190  //Coordinates of the decay positron with respect to the muon spin
191
192  G4double px = stheta*cphi;
193  G4double py = stheta*sphi;
194  G4double pz = ctheta;
195
196  G4ThreeVector direction0(px,py,pz);
197
198  direction0.rotateUz(parent_polarization);
199
200  G4DynamicParticle * daughterparticle0
201    = new G4DynamicParticle( daughters[0], daughtermomentum[0]*direction0);
202
203  products->PushProducts(daughterparticle0);
204
205
206  // daughter 1 ,2 (neutrinos)
207  // create neutrinos in the C.M frame of two neutrinos
208  G4double energy2 = parentmass*(1.0 - x/2.0); 
209  G4double vmass   = std::sqrt((energy2-daughtermomentum[0])*(energy2+daughtermomentum[0]));
210  G4double beta = -1.0*daughtermomentum[0]/energy2;
211  G4double costhetan = 2.*G4UniformRand()-1.0;
212  G4double sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
213  G4double phin  = twopi*G4UniformRand()*rad;
214  G4double sinphin = std::sin(phin);
215  G4double cosphin = std::cos(phin);
216
217  G4ThreeVector direction1(sinthetan*cosphin,sinthetan*sinphin,costhetan);
218  G4DynamicParticle * daughterparticle1
219    = new G4DynamicParticle( daughters[1], direction1*(vmass/2.));
220  G4DynamicParticle * daughterparticle2
221    = new G4DynamicParticle( daughters[2], direction1*(-1.0*vmass/2.));
222
223  // boost to the muon rest frame
224  G4LorentzVector p4;
225  p4 = daughterparticle1->Get4Momentum();
226  p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
227  daughterparticle1->Set4Momentum(p4);
228  p4 = daughterparticle2->Get4Momentum();
229  p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
230  daughterparticle2->Set4Momentum(p4);
231  products->PushProducts(daughterparticle1);
232  products->PushProducts(daughterparticle2);
233  daughtermomentum[1] = daughterparticle1->GetTotalMomentum();
234  daughtermomentum[2] = daughterparticle2->GetTotalMomentum();
235
236  // output message
237#ifdef G4VERBOSE
238  if (GetVerboseLevel()>1) {
239    G4cout << "G4MuonDecayChannelWithSpin::DecayIt ";
240    G4cout << "  create decay products in rest frame " <<G4endl;
241    products->DumpInfo();
242  }
243#endif
244  return products;
245}
246
247G4double G4MuonDecayChannelWithSpin::R_c(G4double x){
248
249  G4int n_max = (int)(100.*x);
250
251  if(n_max<10)n_max=10;
252
253  G4double L2 = 0.0;
254
255  for(G4int n=1; n<=n_max; n++){
256    L2 += std::pow(x,n)/(n*n);
257  }
258
259  G4double omega = std::log(EMMU/EMASS);
260
261  G4double r_c;
262
263  r_c = 2.*L2-(pi*pi/3.)-2.;
264  r_c = r_c + omega * (1.5+2.*std::log((1.-x)/x));
265  r_c = r_c - std::log(x)*(2.*std::log(x)-1.);
266  r_c = r_c + (3.*std::log(x)-1.-1./x)*std::log(1.-x);
267
268  return r_c;
269}
Note: See TracBrowser for help on using the repository browser.