| [824] | 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 | // 25 July 2007 P.Gumplinger
|
|---|
| 31 | // Samples Radiative Muon Decay
|
|---|
| 32 | // Reference:
|
|---|
| 33 | // TRIUMF/TWIST Technote:
|
|---|
| 34 | // "Radiative muon decay" by P. Depommier and A. Vacheret
|
|---|
| 35 | //
|
|---|
| 36 | // ------------------------------------------------------------
|
|---|
| 37 | //
|
|---|
| 38 | //
|
|---|
| 39 | //
|
|---|
| 40 | //
|
|---|
| 41 |
|
|---|
| 42 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 43 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 44 |
|
|---|
| 45 | #ifndef G4MuonRadiativeDecayChannelWithSpin_h
|
|---|
| 46 | #define G4MuonRadiativeDecayChannelWithSpin_h 1
|
|---|
| 47 |
|
|---|
| 48 | #include "globals.hh"
|
|---|
| 49 | #include "Randomize.hh"
|
|---|
| 50 | #include "G4ThreeVector.hh"
|
|---|
| 51 | #include "G4VDecayChannel.hh"
|
|---|
| 52 |
|
|---|
| 53 | class G4MuonRadiativeDecayChannelWithSpin : public G4VDecayChannel
|
|---|
| 54 | {
|
|---|
| 55 | // Class Decription
|
|---|
| 56 | // This class describes radiative muon decay kinemtics, but
|
|---|
| 57 | // gives incorrect energy spectrum for neutrinos
|
|---|
| 58 |
|
|---|
| 59 | public: // With Description
|
|---|
| 60 |
|
|---|
| 61 | //Constructors
|
|---|
| 62 | G4MuonRadiativeDecayChannelWithSpin(const G4String& theParentName,
|
|---|
| 63 | G4double theBR);
|
|---|
| 64 | // Destructor
|
|---|
| 65 | virtual ~G4MuonRadiativeDecayChannelWithSpin();
|
|---|
| 66 |
|
|---|
| 67 | public: // With Description
|
|---|
| 68 |
|
|---|
| 69 | virtual G4DecayProducts *DecayIt(G4double);
|
|---|
| 70 |
|
|---|
| 71 | void SetPolarization(G4ThreeVector);
|
|---|
| 72 | const G4ThreeVector& GetPolarization() const;
|
|---|
| 73 |
|
|---|
| 74 | private:
|
|---|
| 75 |
|
|---|
| 76 | G4ThreeVector parent_polarization;
|
|---|
| 77 |
|
|---|
| 78 | G4double fron(G4double Pmu, G4double x, G4double y,
|
|---|
| 79 | G4double cthetaE, G4double cthetaG, G4double cthetaEG);
|
|---|
| 80 |
|
|---|
| 81 | // rn3dim generates random vectors, uniformly distributed over the surface
|
|---|
| 82 | // of a sphere of given radius.
|
|---|
| 83 | // See http://wwwinfo.cern.ch/asdoc/shortwrupsdir/v131/top.html
|
|---|
| 84 |
|
|---|
| 85 | void rn3dim(G4double& x, G4double& y, G4double& z, G4double xlong);
|
|---|
| 86 |
|
|---|
| 87 | G4double atan4(G4double x, G4double y);
|
|---|
| 88 |
|
|---|
| 89 | G4double EMMU;
|
|---|
| 90 | G4double EMASS;
|
|---|
| 91 |
|
|---|
| 92 | };
|
|---|
| 93 |
|
|---|
| 94 | inline void G4MuonRadiativeDecayChannelWithSpin::
|
|---|
| 95 | SetPolarization(G4ThreeVector polar)
|
|---|
| 96 | {
|
|---|
| 97 | parent_polarization = polar;
|
|---|
| 98 | }
|
|---|
| 99 |
|
|---|
| 100 | inline const G4ThreeVector& G4MuonRadiativeDecayChannelWithSpin::
|
|---|
| 101 | GetPolarization() const
|
|---|
| 102 | {
|
|---|
| 103 | return parent_polarization;
|
|---|
| 104 | }
|
|---|
| 105 |
|
|---|
| 106 | inline void G4MuonRadiativeDecayChannelWithSpin::rn3dim(G4double& x,
|
|---|
| 107 | G4double& y,
|
|---|
| 108 | G4double& z,
|
|---|
| 109 | G4double xlong)
|
|---|
| 110 | {
|
|---|
| 111 | G4double a = 0.; G4double b = 0.; G4double c = 0.; G4double r = 0.;
|
|---|
| 112 |
|
|---|
| 113 | do {
|
|---|
| 114 | a = G4UniformRand() - 0.5;
|
|---|
| 115 | b = G4UniformRand() - 0.5;
|
|---|
| 116 | c = G4UniformRand() - 0.5;
|
|---|
| 117 | r = a*a + b*b + c*c;
|
|---|
| 118 | } while (r > 0.25);
|
|---|
| 119 |
|
|---|
| 120 | G4double rinv = xlong/(std::sqrt(r));
|
|---|
| 121 | x = a * rinv;
|
|---|
| 122 | y = b * rinv;
|
|---|
| 123 | z = c * rinv;
|
|---|
| 124 |
|
|---|
| 125 | return;
|
|---|
| 126 | }
|
|---|
| 127 |
|
|---|
| 128 | inline G4double G4MuonRadiativeDecayChannelWithSpin::atan4(G4double x,
|
|---|
| 129 | G4double y)
|
|---|
| 130 | {
|
|---|
| 131 | G4double phi = 0.;
|
|---|
| 132 |
|
|---|
| 133 | G4double pi = 4.*std::atan(1.);
|
|---|
| 134 |
|
|---|
| 135 | if (x==0. && y>0.){
|
|---|
| 136 | phi = 0.5*pi;
|
|---|
| 137 | } else if (x==0. && y<0.){
|
|---|
| 138 | phi = 1.5*pi;
|
|---|
| 139 | } else if (y==0. && x>0.){
|
|---|
| 140 | phi = 0.;
|
|---|
| 141 | } else if (y==0. && x<0.){
|
|---|
| 142 | phi = pi;
|
|---|
| 143 | } else if (x>0. ){
|
|---|
| 144 | phi = std::atan(y/x);
|
|---|
| 145 | } else if (x<0. ){
|
|---|
| 146 | phi = std::atan(y/x) + pi;
|
|---|
| 147 | }
|
|---|
| 148 |
|
|---|
| 149 | return phi;
|
|---|
| 150 | }
|
|---|
| 151 |
|
|---|
| 152 | #endif
|
|---|