Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuonDecayChannelWithSpin Class Reference

#include <G4MuonDecayChannelWithSpin.hh>

+ Inheritance diagram for G4MuonDecayChannelWithSpin:

Public Member Functions

 G4MuonDecayChannelWithSpin (const G4String &theParentName, G4double theBR)
 
virtual ~G4MuonDecayChannelWithSpin ()
 
virtual G4DecayProductsDecayIt (G4double)
 
void SetPolarization (G4ThreeVector)
 
const G4ThreeVectorGetPolarization () const
 
- Public Member Functions inherited from G4MuonDecayChannel
 G4MuonDecayChannel (const G4String &theParentName, G4double theBR)
 
virtual ~G4MuonDecayChannel ()
 
virtual G4DecayProductsDecayIt (G4double)
 
- Public Member Functions inherited from G4VDecayChannel
 G4VDecayChannel (const G4String &aName, G4int Verbose=1)
 
 G4VDecayChannel (const G4String &aName, const G4String &theParentName, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="", const G4String &theDaughterName4="")
 
virtual ~G4VDecayChannel ()
 
G4int operator== (const G4VDecayChannel &right) const
 
G4int operator!= (const G4VDecayChannel &right) const
 
G4int operator< (const G4VDecayChannel &right) const
 
virtual G4DecayProductsDecayIt (G4double parentMass=-1.0)=0
 
const G4StringGetKinematicsName () const
 
G4double GetBR () const
 
G4int GetNumberOfDaughters () const
 
G4ParticleDefinitionGetParent ()
 
G4ParticleDefinitionGetDaughter (G4int anIndex)
 
G4int GetAngularMomentum ()
 
const G4StringGetParentName () const
 
const G4StringGetDaughterName (G4int anIndex) const
 
G4double GetParentMass () const
 
G4double GetDaughterMass (G4int anIndex) const
 
void SetParent (const G4ParticleDefinition *particle_type)
 
void SetParent (const G4String &particle_name)
 
void SetBR (G4double value)
 
void SetNumberOfDaughters (G4int value)
 
void SetDaughter (G4int anIndex, const G4ParticleDefinition *particle_type)
 
void SetDaughter (G4int anIndex, const G4String &particle_name)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void DumpInfo ()
 

Protected Member Functions

 G4MuonDecayChannelWithSpin (const G4MuonDecayChannelWithSpin &)
 
G4MuonDecayChannelWithSpinoperator= (const G4MuonDecayChannelWithSpin &)
 
- Protected Member Functions inherited from G4MuonDecayChannel
 G4MuonDecayChannel (const G4MuonDecayChannel &)
 
G4MuonDecayChanneloperator= (const G4MuonDecayChannel &)
 
 G4MuonDecayChannel ()
 
- Protected Member Functions inherited from G4VDecayChannel
 G4VDecayChannel ()
 
 G4VDecayChannel (const G4VDecayChannel &)
 
G4VDecayChanneloperator= (const G4VDecayChannel &)
 
void ClearDaughtersName ()
 
void FillDaughters ()
 
void FillParent ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VDecayChannel
G4String kinematics_name
 
G4double rbranch
 
G4int numberOfDaughters
 
G4Stringparent_name
 
G4String ** daughters_name
 
G4ParticleTableparticletable
 
G4ParticleDefinitionparent
 
G4ParticleDefinition ** daughters
 
G4double parent_mass
 
G4doubledaughters_mass
 
G4int verboseLevel
 
- Static Protected Attributes inherited from G4VDecayChannel
static const G4String noName = " "
 

Detailed Description

Definition at line 50 of file G4MuonDecayChannelWithSpin.hh.

Constructor & Destructor Documentation

◆ G4MuonDecayChannelWithSpin() [1/2]

G4MuonDecayChannelWithSpin::G4MuonDecayChannelWithSpin ( const G4String theParentName,
G4double  theBR 
)

Definition at line 59 of file G4MuonDecayChannelWithSpin.cc.

61 : G4MuonDecayChannel(theParentName,theBR),
62 parent_polarization(),
63 EMMU( 0.*MeV),
64 EMASS( 0.*MeV)
65{
66}

Referenced by G4MuonDecayChannelWithSpin().

◆ ~G4MuonDecayChannelWithSpin()

G4MuonDecayChannelWithSpin::~G4MuonDecayChannelWithSpin ( )
virtual

Definition at line 68 of file G4MuonDecayChannelWithSpin.cc.

69{
70}

◆ G4MuonDecayChannelWithSpin() [2/2]

G4MuonDecayChannelWithSpin::G4MuonDecayChannelWithSpin ( const G4MuonDecayChannelWithSpin right)
protected

Definition at line 72 of file G4MuonDecayChannelWithSpin.cc.

72 :
74{
75 parent_polarization = right.parent_polarization;
76 EMMU = right.EMMU;
77 EMASS = right.EMASS;
78}

Member Function Documentation

◆ DecayIt()

G4DecayProducts * G4MuonDecayChannelWithSpin::DecayIt ( G4double  )
virtual

Reimplemented from G4MuonDecayChannel.

Definition at line 111 of file G4MuonDecayChannelWithSpin.cc.

112{
113 // This version assumes V-A coupling with 1st order radiative correctons,
114 // the standard model Michel parameter values, but
115 // gives incorrect energy spectrum for neutrinos
116
117#ifdef G4VERBOSE
118 if (GetVerboseLevel()>1) G4cout << "G4MuonDecayChannelWithSpin::DecayIt ";
119#endif
120
121 if (parent == 0) FillParent();
122 if (daughters == 0) FillDaughters();
123
124 // parent mass
125 G4double parentmass = parent->GetPDGMass();
126
127 EMMU = parentmass;
128
129 //daughters'mass
130 G4double daughtermass[3];
131 G4double sumofdaughtermass = 0.0;
132 for (G4int index=0; index<3; index++){
133 daughtermass[index] = daughters[index]->GetPDGMass();
134 sumofdaughtermass += daughtermass[index];
135 }
136
137 EMASS = daughtermass[0];
138
139 //create parent G4DynamicParticle at rest
140 G4ThreeVector dummy;
141 G4DynamicParticle * parentparticle = new G4DynamicParticle( parent, dummy, 0.0);
142 //create G4Decayproducts
143 G4DecayProducts *products = new G4DecayProducts(*parentparticle);
144 delete parentparticle;
145
146 // calcurate electron energy
147
148 G4double michel_rho = 0.75; //Standard Model Michel rho
149 G4double michel_delta = 0.75; //Standard Model Michel delta
150 G4double michel_xsi = 1.00; //Standard Model Michel xsi
151 G4double michel_eta = 0.00; //Standard Model eta
152
153 G4double rndm, x, ctheta;
154
155 G4double FG;
156 G4double FG_max = 2.00;
157
158 G4double W_mue = (EMMU*EMMU+EMASS*EMASS)/(2.*EMMU);
159 G4double x0 = EMASS/W_mue;
160
161 G4double x0_squared = x0*x0;
162
163 // ***************************************************
164 // x0 <= x <= 1. and -1 <= y <= 1
165 //
166 // F(x,y) = f(x)*g(x,y); g(x,y) = 1.+g(x)*y
167 // ***************************************************
168
169 // ***** sampling F(x,y) directly (brute force) *****
170
171 do{
172
173 // Sample the positron energy by sampling from F
174
175 rndm = G4UniformRand();
176
177 x = x0 + rndm*(1.-x0);
178
179 G4double x_squared = x*x;
180
181 G4double F_IS, F_AS, G_IS, G_AS;
182
183 F_IS = 1./6.*(-2.*x_squared+3.*x-x0_squared);
184 F_AS = 1./6.*std::sqrt(x_squared-x0_squared)*(2.*x-2.+std::sqrt(1.-x0_squared));
185
186 G_IS = 2./9.*(michel_rho-0.75)*(4.*x_squared-3.*x-x0_squared);
187 G_IS = G_IS + michel_eta*(1.-x)*x0;
188
189 G_AS = 3.*(michel_xsi-1.)*(1.-x);
190 G_AS = G_AS+2.*(michel_xsi*michel_delta-0.75)*(4.*x-4.+std::sqrt(1.-x0_squared));
191 G_AS = 1./9.*std::sqrt(x_squared-x0_squared)*G_AS;
192
193 F_IS = F_IS + G_IS;
194 F_AS = F_AS + G_AS;
195
196 // *** Radiative Corrections ***
197
198 G4double R_IS = F_c(x,x0);
199
200 G4double F = 6.*F_IS + R_IS/std::sqrt(x_squared-x0_squared);
201
202 // *** Radiative Corrections ***
203
204 G4double R_AS = F_theta(x,x0);
205
206 rndm = G4UniformRand();
207
208 ctheta = 2.*rndm-1.;
209
210 G4double G = 6.*F_AS - R_AS/std::sqrt(x_squared-x0_squared);
211
212 FG = std::sqrt(x_squared-x0_squared)*F*(1.+(G/F)*ctheta);
213
214 if(FG>FG_max){
215 G4cout<<"***Problem in Muon Decay *** : FG > FG_max"<<G4endl;
216 FG_max = FG;
217 }
218
219 rndm = G4UniformRand();
220
221 }while(FG<rndm*FG_max);
222
223 G4double energy = x * W_mue;
224
225 rndm = G4UniformRand();
226
227 G4double phi = twopi * rndm;
228
229 if(energy < EMASS) energy = EMASS;
230
231 // calculate daughter momentum
232 G4double daughtermomentum[3];
233
234 daughtermomentum[0] = std::sqrt(energy*energy - EMASS*EMASS);
235
236 G4double stheta = std::sqrt(1.-ctheta*ctheta);
237 G4double cphi = std::cos(phi);
238 G4double sphi = std::sin(phi);
239
240 //Coordinates of the decay positron with respect to the muon spin
241
242 G4double px = stheta*cphi;
243 G4double py = stheta*sphi;
244 G4double pz = ctheta;
245
246 G4ThreeVector direction0(px,py,pz);
247
248 direction0.rotateUz(parent_polarization);
249
250 G4DynamicParticle * daughterparticle0
251 = new G4DynamicParticle( daughters[0], daughtermomentum[0]*direction0);
252
253 products->PushProducts(daughterparticle0);
254
255
256 // daughter 1 ,2 (neutrinos)
257 // create neutrinos in the C.M frame of two neutrinos
258 G4double energy2 = parentmass*(1.0 - x/2.0);
259 G4double vmass = std::sqrt((energy2-daughtermomentum[0])*(energy2+daughtermomentum[0]));
260 G4double beta = -1.0*daughtermomentum[0]/energy2;
261 G4double costhetan = 2.*G4UniformRand()-1.0;
262 G4double sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
263 G4double phin = twopi*G4UniformRand()*rad;
264 G4double sinphin = std::sin(phin);
265 G4double cosphin = std::cos(phin);
266
267 G4ThreeVector direction1(sinthetan*cosphin,sinthetan*sinphin,costhetan);
268 G4DynamicParticle * daughterparticle1
269 = new G4DynamicParticle( daughters[1], direction1*(vmass/2.));
270 G4DynamicParticle * daughterparticle2
271 = new G4DynamicParticle( daughters[2], direction1*(-1.0*vmass/2.));
272
273 // boost to the muon rest frame
275 p4 = daughterparticle1->Get4Momentum();
276 p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
277 daughterparticle1->Set4Momentum(p4);
278 p4 = daughterparticle2->Get4Momentum();
279 p4.boost( direction0.x()*beta, direction0.y()*beta, direction0.z()*beta);
280 daughterparticle2->Set4Momentum(p4);
281 products->PushProducts(daughterparticle1);
282 products->PushProducts(daughterparticle2);
283 daughtermomentum[1] = daughterparticle1->GetTotalMomentum();
284 daughtermomentum[2] = daughterparticle2->GetTotalMomentum();
285
286 // output message
287#ifdef G4VERBOSE
288 if (GetVerboseLevel()>1) {
289 G4cout << "G4MuonDecayChannelWithSpin::DecayIt ";
290 G4cout << " create decay products in rest frame " <<G4endl;
291 products->DumpInfo();
292 }
293#endif
294 return products;
295}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
HepLorentzVector & boost(double, double, double)
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetTotalMomentum() const
G4int GetVerboseLevel() const
G4ParticleDefinition * parent
G4ParticleDefinition ** daughters

◆ GetPolarization()

const G4ThreeVector & G4MuonDecayChannelWithSpin::GetPolarization ( ) const
inline

Definition at line 101 of file G4MuonDecayChannelWithSpin.hh.

102{
103 return parent_polarization;
104}

◆ operator=()

G4MuonDecayChannelWithSpin & G4MuonDecayChannelWithSpin::operator= ( const G4MuonDecayChannelWithSpin right)
protected

Definition at line 80 of file G4MuonDecayChannelWithSpin.cc.

81{
82 if (this != &right) {
85 rbranch = right.rbranch;
86
87 // copy parent name
88 parent_name = new G4String(*right.parent_name);
89
90 // clear daughters_name array
92
93 // recreate array
95 if ( numberOfDaughters >0 ) {
98 //copy daughters name
99 for (G4int index=0; index < numberOfDaughters; index++) {
100 daughters_name[index] = new G4String(*right.daughters_name[index]);
101 }
102 }
103 parent_polarization = right.parent_polarization;
104 EMMU = right.EMMU;
105 EMASS = right.EMASS;
106 }
107 return *this;
108}
G4String * parent_name
G4String ** daughters_name
G4String kinematics_name

◆ SetPolarization()

void G4MuonDecayChannelWithSpin::SetPolarization ( G4ThreeVector  polar)
inline

Definition at line 96 of file G4MuonDecayChannelWithSpin.hh.

97{
98 parent_polarization = polar;
99}

Referenced by G4DecayWithSpin::DecayIt().


The documentation for this class was generated from the following files: