Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TauLeptonicDecayChannel.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// G4TauLeptonicDecayChannel class implementation
27//
28// Author: H.Kurashige, 30 May 1997
29// --------------------------------------------------------------------
30
33#include "G4SystemOfUnits.hh"
34#include "G4DecayProducts.hh"
35#include "G4VDecayChannel.hh"
37#include "Randomize.hh"
38#include "G4LorentzVector.hh"
39#include "G4LorentzRotation.hh"
40
41
44{
45}
46
47// --------------------------------------------------------------------
49G4TauLeptonicDecayChannel(const G4String& theParentName,
50 G4double theBR,
51 const G4String& theLeptonName)
52 : G4VDecayChannel("Tau Leptonic Decay", 1)
53{
54 // set names for daughter particles
55 if (theParentName == "tau+")
56 {
57 SetBR(theBR);
58 SetParent("tau+");
60 if ((theLeptonName=="e-"||theLeptonName=="e+"))
61 {
62 SetDaughter(0, "e+");
63 SetDaughter(1, "nu_e");
64 SetDaughter(2, "anti_nu_tau");
65 }
66 else
67 {
68 SetDaughter(0, "mu+");
69 SetDaughter(1, "nu_mu");
70 SetDaughter(2, "anti_nu_tau");
71 }
72 }
73 else if (theParentName == "tau-")
74 {
75 SetBR(theBR);
76 SetParent("tau-");
78 if ((theLeptonName=="e-"||theLeptonName=="e+"))
79 {
80 SetDaughter(0, "e-");
81 SetDaughter(1, "anti_nu_e");
82 SetDaughter(2, "nu_tau");
83 }
84 else
85 {
86 SetDaughter(0, "mu-");
87 SetDaughter(1, "anti_nu_mu");
88 SetDaughter(2, "nu_tau");
89 }
90 }
91 else
92 {
93#ifdef G4VERBOSE
94 if (GetVerboseLevel()>0)
95 {
96 G4cout << "G4TauLeptonicDecayChannel:: constructor :";
97 G4cout << " parent particle is not tau but ";
98 G4cout << theParentName << G4endl;
99 }
100#endif
101 }
102}
103
104// --------------------------------------------------------------------
106{
107}
108
109// --------------------------------------------------------------------
112 : G4VDecayChannel(right)
113{
114}
115
116// --------------------------------------------------------------------
119{
120 if (this != &right)
121 {
124 rbranch = right.rbranch;
125
126 // copy parent name
127 parent_name = new G4String(*right.parent_name);
128
129 // clear daughters_name array
131
132 // recreate array
134 if ( numberOfDaughters >0 )
135 {
138 // copy daughters name
139 for (G4int index=0; index<numberOfDaughters; ++index)
140 {
141 daughters_name[index] = new G4String(*right.daughters_name[index]);
142 }
143 }
144 }
145 return *this;
146}
147
148// --------------------------------------------------------------------
150{
151 // this version neglects muon polarization
152 // assumes the pure V-A coupling
153 // gives incorrect energy spectrum for neutrinos
154
155#ifdef G4VERBOSE
156 if (GetVerboseLevel()>1) G4cout << "G4TauLeptonicDecayChannel::DecayIt()";
157#endif
158
161
162 // parent mass
163 G4double parentmass = G4MT_parent->GetPDGMass();
164
165 // daughters'mass
166 const G4int N_DAUGHTER=3;
167 G4double daughtermass[N_DAUGHTER];
168 for (G4int index=0; index<N_DAUGHTER; ++index)
169 {
170 daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
171 }
172
173 // create parent G4DynamicParticle at rest
174 G4ThreeVector dummy;
175 G4DynamicParticle* parentparticle
176 = new G4DynamicParticle(G4MT_parent, dummy, 0.0);
177 // create G4Decayproducts
178 G4DecayProducts* products = new G4DecayProducts(*parentparticle);
179 delete parentparticle;
180
181 // calculate daughter momentum
182 G4double daughtermomentum[N_DAUGHTER];
183
184 // calculate lepton momentum
185 G4double pmax = (parentmass*parentmass
186 - daughtermass[0]*daughtermass[0])/2./parentmass;
187 G4double p, e;
188 G4double r;
189 const std::size_t MAX_LOOP=10000;
190 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
191 {
192 // determine momentum/energy
193 r = G4UniformRand();
194 p = pmax*G4UniformRand();
195 e = std::sqrt(p*p + daughtermass[0]*daughtermass[0]);
196 if (r < spectrum(p,e,parentmass,daughtermass[0]) ) break;
197 }
198
199 // create daughter G4DynamicParticle
200 // daughter 0 (lepton)
201 daughtermomentum[0] = p;
202 G4double costheta, sintheta, phi, sinphi, cosphi;
203 costheta = 2.*G4UniformRand()-1.0;
204 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
205 phi = twopi*G4UniformRand()*rad;
206 sinphi = std::sin(phi);
207 cosphi = std::cos(phi);
208 G4ThreeVector direction0(sintheta*cosphi,sintheta*sinphi,costheta);
209 G4DynamicParticle * daughterparticle
210 = new G4DynamicParticle(G4MT_daughters[0], direction0*daughtermomentum[0]);
211 products->PushProducts(daughterparticle);
212
213 // daughter 1 ,2 (nutrinos)
214 // create neutrinos in the C.M frame of two neutrinos
215 G4double energy2 = parentmass-e;
216 G4double vmass = std::sqrt((energy2-daughtermomentum[0])
217 *(energy2+daughtermomentum[0]));
218 G4double beta = -1.0*daughtermomentum[0]/energy2;
219 G4double costhetan = 2.*G4UniformRand()-1.0;
220 G4double sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
221 G4double phin = twopi*G4UniformRand()*rad;
222 G4double sinphin = std::sin(phin);
223 G4double cosphin = std::cos(phin);
224
225 G4ThreeVector direction1(sinthetan*cosphin,sinthetan*sinphin,costhetan);
226 G4DynamicParticle * daughterparticle1
227 = new G4DynamicParticle( G4MT_daughters[1], direction1*(vmass/2.));
228 G4DynamicParticle * daughterparticle2
229 = new G4DynamicParticle( G4MT_daughters[2], direction1*(-1.0*vmass/2.));
230
231 // boost to the muon rest frame
233 p4 = daughterparticle1->Get4Momentum();
234 p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
235 daughterparticle1->Set4Momentum(p4);
236 p4 = daughterparticle2->Get4Momentum();
237 p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
238 daughterparticle2->Set4Momentum(p4);
239 products->PushProducts(daughterparticle1);
240 products->PushProducts(daughterparticle2);
241 daughtermomentum[1] = daughterparticle1->GetTotalMomentum();
242 daughtermomentum[2] = daughterparticle2->GetTotalMomentum();
243
244 // output message
245#ifdef G4VERBOSE
246 if (GetVerboseLevel()>1)
247 {
248 G4cout << "G4TauLeptonicDecayChannel::DecayIt ";
249 G4cout << " create decay products in rest frame " << G4endl;
250 products->DumpInfo();
251 }
252#endif
253 return products;
254}
255
256// --------------------------------------------------------------------
257G4double G4TauLeptonicDecayChannel::spectrum(G4double p,
258 G4double e,
259 G4double mtau,
260 G4double ml)
261{
262 G4double f1;
263 f1 = 3.0*e*(mtau*mtau+ml*ml)-4.0*mtau*e*e-2.0*mtau*ml*ml;
264 return p*(f1)/(mtau*mtau*mtau*mtau)/(0.6);
265}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
HepLorentzVector & boost(double, double, double)
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetTotalMomentum() const
G4TauLeptonicDecayChannel & operator=(const G4TauLeptonicDecayChannel &)
virtual G4DecayProducts * DecayIt(G4double)
G4ParticleDefinition ** G4MT_daughters
G4String * parent_name
void SetBR(G4double value)
G4String ** daughters_name
G4int GetVerboseLevel() const
void SetNumberOfDaughters(G4int value)
G4ParticleDefinition * G4MT_parent
void CheckAndFillDaughters()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
G4String kinematics_name
void SetParent(const G4ParticleDefinition *particle_type)