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 | #ifndef G4MuonDecayChannelWithSpin_hh |
---|
42 | #define G4MuonDecayChannelWithSpin_hh 1 |
---|
43 | |
---|
44 | #include "globals.hh" |
---|
45 | #include "G4ThreeVector.hh" |
---|
46 | #include "G4MuonDecayChannel.hh" |
---|
47 | |
---|
48 | class G4MuonDecayChannelWithSpin : public G4MuonDecayChannel |
---|
49 | { |
---|
50 | // Class Decription |
---|
51 | // This class describes muon decay kinemtics. |
---|
52 | // This version assumes V-A coupling with 1st order radiative correctons, |
---|
53 | // the standard model Michel parameter values, but |
---|
54 | // gives incorrect energy spectrum for neutrinos |
---|
55 | |
---|
56 | public: // With Description |
---|
57 | |
---|
58 | //Constructors |
---|
59 | G4MuonDecayChannelWithSpin(const G4String& theParentName, |
---|
60 | G4double theBR); |
---|
61 | // Destructor |
---|
62 | virtual ~G4MuonDecayChannelWithSpin(); |
---|
63 | |
---|
64 | public: // With Description |
---|
65 | |
---|
66 | virtual G4DecayProducts *DecayIt(G4double); |
---|
67 | |
---|
68 | void SetPolarization(G4ThreeVector); |
---|
69 | const G4ThreeVector& GetPolarization() const; |
---|
70 | |
---|
71 | private: |
---|
72 | |
---|
73 | G4ThreeVector parent_polarization; |
---|
74 | |
---|
75 | // Radiative Correction Factors |
---|
76 | |
---|
77 | G4double F_c(G4double x, G4double x0); |
---|
78 | G4double F_theta(G4double x, G4double x0); |
---|
79 | G4double R_c(G4double x); |
---|
80 | |
---|
81 | G4double EMMU; |
---|
82 | G4double EMASS; |
---|
83 | |
---|
84 | }; |
---|
85 | |
---|
86 | inline void G4MuonDecayChannelWithSpin::SetPolarization(G4ThreeVector polar) |
---|
87 | { |
---|
88 | parent_polarization = polar; |
---|
89 | } |
---|
90 | |
---|
91 | inline const G4ThreeVector& G4MuonDecayChannelWithSpin::GetPolarization() const |
---|
92 | { |
---|
93 | return parent_polarization; |
---|
94 | } |
---|
95 | |
---|
96 | inline G4double G4MuonDecayChannelWithSpin::F_c(G4double x, G4double x0) |
---|
97 | { |
---|
98 | G4double omega = std::log(EMMU/EMASS); |
---|
99 | |
---|
100 | G4double f_c; |
---|
101 | |
---|
102 | f_c = (5.+17.*x-34.*x*x)*(omega+std::log(x))-22.*x+34.*x*x; |
---|
103 | f_c = (1.-x)/(3.*x*x)*f_c; |
---|
104 | f_c = (6.-4.*x)*R_c(x)+(6.-6.*x)*std::log(x) + f_c; |
---|
105 | f_c = (fine_structure_const/twopi) * (x*x-x0*x0) * f_c; |
---|
106 | |
---|
107 | return f_c; |
---|
108 | } |
---|
109 | |
---|
110 | inline G4double G4MuonDecayChannelWithSpin::F_theta(G4double x, G4double x0) |
---|
111 | { |
---|
112 | G4double omega = std::log(EMMU/EMASS); |
---|
113 | |
---|
114 | G4double f_theta; |
---|
115 | |
---|
116 | f_theta = (1.+x+34*x*x)*(omega+std::log(x))+3.-7.*x-32.*x*x; |
---|
117 | f_theta = f_theta + ((4.*(1.-x)*(1.-x))/x)*std::log(1.-x); |
---|
118 | f_theta = (1.-x)/(3.*x*x) * f_theta; |
---|
119 | f_theta = (2.-4.*x)*R_c(x)+(2.-6.*x)*std::log(x)-f_theta; |
---|
120 | f_theta = (fine_structure_const/twopi) * (x*x-x0*x0) * f_theta; |
---|
121 | |
---|
122 | return f_theta; |
---|
123 | } |
---|
124 | |
---|
125 | #endif |
---|