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