Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions
G4MuonDecayChannelWithSpin Class Reference

#include <G4MuonDecayChannelWithSpin.hh>

Inheritance diagram for G4MuonDecayChannelWithSpin:
G4MuonDecayChannel G4VDecayChannel

Public Member Functions

virtual G4DecayProductsDecayIt (G4double)
 
void DumpInfo ()
 
 G4MuonDecayChannelWithSpin (const G4String &theParentName, G4double theBR)
 
G4int GetAngularMomentum ()
 
G4double GetBR () const
 
G4ParticleDefinitionGetDaughter (G4int anIndex)
 
G4double GetDaughterMass (G4int anIndex) const
 
const G4StringGetDaughterName (G4int anIndex) const
 
const G4StringGetKinematicsName () const
 
G4int GetNumberOfDaughters () const
 
G4ParticleDefinitionGetParent ()
 
G4double GetParentMass () const
 
const G4StringGetParentName () const
 
const G4ThreeVectorGetPolarization () const
 
G4double GetRangeMass () const
 
G4int GetVerboseLevel () const
 
virtual G4bool IsOKWithParentMass (G4double parentMass)
 
G4bool operator!= (const G4VDecayChannel &r) const
 
G4bool operator< (const G4VDecayChannel &right) const
 
G4bool operator== (const G4VDecayChannel &r) const
 
void SetBR (G4double value)
 
void SetDaughter (G4int anIndex, const G4ParticleDefinition *particle_type)
 
void SetDaughter (G4int anIndex, const G4String &particle_name)
 
void SetNumberOfDaughters (G4int value)
 
void SetParent (const G4ParticleDefinition *particle_type)
 
void SetParent (const G4String &particle_name)
 
void SetPolarization (const G4ThreeVector &)
 
void SetRangeMass (G4double val)
 
void SetVerboseLevel (G4int value)
 
virtual ~G4MuonDecayChannelWithSpin ()
 

Protected Member Functions

void CheckAndFillDaughters ()
 
void CheckAndFillParent ()
 
void ClearDaughtersName ()
 
G4double DynamicalMass (G4double massPDG, G4double width, G4double maxDev=1.0) const
 
 G4MuonDecayChannelWithSpin (const G4MuonDecayChannelWithSpin &)
 
G4MuonDecayChannelWithSpinoperator= (const G4MuonDecayChannelWithSpin &)
 

Protected Attributes

G4String ** daughters_name = nullptr
 
G4Mutex daughtersMutex
 
G4ParticleDefinition ** G4MT_daughters = nullptr
 
G4doubleG4MT_daughters_mass = nullptr
 
G4doubleG4MT_daughters_width = nullptr
 
G4ParticleDefinitionG4MT_parent = nullptr
 
G4double G4MT_parent_mass = 0.0
 
G4String kinematics_name = ""
 
G4int numberOfDaughters = 0
 
G4Stringparent_name = nullptr
 
G4ThreeVector parent_polarization
 
G4Mutex parentMutex
 
G4ParticleTableparticletable = nullptr
 
G4double rangeMass = 2.5
 
G4double rbranch = 0.0
 
G4int verboseLevel = 1
 

Static Protected Attributes

static const G4String noName = " "
 

Private Member Functions

G4double F_c (G4double x, G4double x0, G4double omega)
 
G4double F_theta (G4double x, G4double x0, G4double omega)
 
void FillDaughters ()
 
void FillParent ()
 
 G4MuonDecayChannelWithSpin ()
 
const G4StringGetNoName () const
 
G4double R_c (G4double x, G4double omega)
 

Detailed Description

Definition at line 51 of file G4MuonDecayChannelWithSpin.hh.

Constructor & Destructor Documentation

◆ G4MuonDecayChannelWithSpin() [1/3]

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

Definition at line 51 of file G4MuonDecayChannelWithSpin.cc.

54 : G4MuonDecayChannel(theParentName,theBR)
55{
56}

◆ ~G4MuonDecayChannelWithSpin()

G4MuonDecayChannelWithSpin::~G4MuonDecayChannelWithSpin ( )
virtual

Definition at line 58 of file G4MuonDecayChannelWithSpin.cc.

59{
60}

◆ G4MuonDecayChannelWithSpin() [2/3]

G4MuonDecayChannelWithSpin::G4MuonDecayChannelWithSpin ( const G4MuonDecayChannelWithSpin right)
protected

Definition at line 62 of file G4MuonDecayChannelWithSpin.cc.

64 : G4MuonDecayChannel(right)
65{
66}

◆ G4MuonDecayChannelWithSpin() [3/3]

G4MuonDecayChannelWithSpin::G4MuonDecayChannelWithSpin ( )
private

Definition at line 46 of file G4MuonDecayChannelWithSpin.cc.

48{
49}

Member Function Documentation

◆ CheckAndFillDaughters()

void G4VDecayChannel::CheckAndFillDaughters ( )
inlineprotectedinherited

◆ CheckAndFillParent()

void G4VDecayChannel::CheckAndFillParent ( )
inlineprotectedinherited

◆ ClearDaughtersName()

void G4VDecayChannel::ClearDaughtersName ( )
protectedinherited

Definition at line 183 of file G4VDecayChannel.cc.

184{
186 if ( daughters_name != nullptr)
187 {
188 if (numberOfDaughters>0)
189 {
190#ifdef G4VERBOSE
191 if (verboseLevel>1)
192 {
193 G4cout << "G4VDecayChannel::ClearDaughtersName() "
194 << " for " << *parent_name << G4endl;
195 }
196#endif
197 for (G4int index=0; index<numberOfDaughters; ++index)
198 {
199 delete daughters_name[index];
200 }
201 }
202 delete [] daughters_name;
203 daughters_name = nullptr;
204 }
205
206 delete [] G4MT_daughters;
207 delete [] G4MT_daughters_mass;
208 delete [] G4MT_daughters_width;
209 G4MT_daughters_width = nullptr;
210 G4MT_daughters = nullptr;
211 G4MT_daughters_mass = nullptr;
212
214}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4String * parent_name
G4String ** daughters_name
G4double * G4MT_daughters_mass
G4double * G4MT_daughters_width

References G4VDecayChannel::daughters_name, G4VDecayChannel::daughtersMutex, G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, and G4VDecayChannel::verboseLevel.

Referenced by G4DalitzDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4MuonDecayChannel::operator=(), operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4NeutronBetaDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4TauLeptonicDecayChannel::operator=(), G4VDecayChannel::operator=(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::~G4VDecayChannel().

◆ DecayIt()

G4DecayProducts * G4MuonDecayChannelWithSpin::DecayIt ( G4double  )
virtual

Reimplemented from G4MuonDecayChannel.

Definition at line 99 of file G4MuonDecayChannelWithSpin.cc.

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

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, CLHEP::HepLorentzVector::boost(), G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), G4DecayProducts::DumpInfo(), G4INCL::KinematicsUtils::energy(), F_c(), F_theta(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4UniformRand, G4DynamicParticle::Get4Momentum(), G4ParticleDefinition::GetPDGMass(), G4DynamicParticle::GetTotalEnergy(), G4DynamicParticle::GetTotalMomentum(), G4VDecayChannel::GetVerboseLevel(), JustWarning, keV, MeV, G4VDecayChannel::parent_polarization, G4DecayProducts::PushProducts(), rad, CLHEP::Hep3Vector::rotateUz(), G4DynamicParticle::Set4Momentum(), twopi, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DumpInfo()

void G4VDecayChannel::DumpInfo ( )
inherited

Definition at line 560 of file G4VDecayChannel.cc.

561{
562 G4cout << " BR: " << rbranch << " [" << kinematics_name << "]";
563 G4cout << " : " ;
564 for (G4int index=0; index<numberOfDaughters; ++index)
565 {
566 if(daughters_name[index] != nullptr)
567 {
568 G4cout << " " << *(daughters_name[index]);
569 }
570 else
571 {
572 G4cout << " not defined ";
573 }
574 }
575 G4cout << G4endl;
576}
G4String kinematics_name

References G4VDecayChannel::daughters_name, G4cout, G4endl, G4VDecayChannel::kinematics_name, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::rbranch.

Referenced by G4GeneralPhaseSpaceDecay::DecayIt(), G4PhaseSpaceDecayChannel::DecayIt(), and G4KL3DecayChannel::G4KL3DecayChannel().

◆ DynamicalMass()

G4double G4VDecayChannel::DynamicalMass ( G4double  massPDG,
G4double  width,
G4double  maxDev = 1.0 
) const
protectedinherited

Definition at line 585 of file G4VDecayChannel.cc.

587{
588 if (width<=0.0) return massPDG;
589 if (maxDev >rangeMass) maxDev = rangeMass;
590 if (maxDev <=-1.*rangeMass) return massPDG; // cannot calculate
591
592 G4double x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
594 const std::size_t MAX_LOOP=10000;
595 for (std::size_t loop_counter=0; loop_counter <MAX_LOOP; ++loop_counter)
596 {
597 if ( y * (width*width*x*x + massPDG*massPDG*width*width)
598 <= massPDG*massPDG*width*width ) break;
599 x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
600 y = G4UniformRand();
601 }
602 G4double mass = massPDG + x*width;
603 return mass;
604}

References G4UniformRand, and G4VDecayChannel::rangeMass.

Referenced by G4PhaseSpaceDecayChannel::ThreeBodyDecayIt(), and G4PhaseSpaceDecayChannel::TwoBodyDecayIt().

◆ F_c()

G4double G4MuonDecayChannelWithSpin::F_c ( G4double  x,
G4double  x0,
G4double  omega 
)
inlineprivate

Definition at line 83 of file G4MuonDecayChannelWithSpin.hh.

84{
85 G4double f_c;
86
87 f_c = (5.+17.*x-34.*x*x)*(omega+std::log(x))-22.*x+34.*x*x;
88 f_c = (1.-x)/(3.*x*x)*f_c;
89 f_c = (6.-4.*x)*R_c(x,omega)+(6.-6.*x)*std::log(x) + f_c;
90 f_c = (CLHEP::fine_structure_const/CLHEP::twopi) * (x*x-x0*x0) * f_c;
91
92 return f_c;
93}
G4double R_c(G4double x, G4double omega)
static constexpr double fine_structure_const
static constexpr double twopi
Definition: SystemOfUnits.h:56

References CLHEP::fine_structure_const, R_c(), and CLHEP::twopi.

Referenced by DecayIt().

◆ F_theta()

G4double G4MuonDecayChannelWithSpin::F_theta ( G4double  x,
G4double  x0,
G4double  omega 
)
inlineprivate

Definition at line 96 of file G4MuonDecayChannelWithSpin.hh.

97{
98 G4double f_theta;
99
100 f_theta = (1.+x+34*x*x)*(omega+std::log(x))+3.-7.*x-32.*x*x;
101 f_theta = f_theta + ((4.*(1.-x)*(1.-x))/x)*std::log(1.-x);
102 f_theta = (1.-x)/(3.*x*x) * f_theta;
103 f_theta = (2.-4.*x)*R_c(x,omega)+(2.-6.*x)*std::log(x)-f_theta;
104 f_theta = (CLHEP::fine_structure_const/CLHEP::twopi) * (x*x-x0*x0) * f_theta;
105
106 return f_theta;
107}

References CLHEP::fine_structure_const, R_c(), and CLHEP::twopi.

Referenced by DecayIt().

◆ FillDaughters()

void G4VDecayChannel::FillDaughters ( )
privateinherited

Definition at line 308 of file G4VDecayChannel.cc.

309{
311
312 // Double check, check again if another thread has already filled this, in
313 // case do not need to do anything
314 if ( G4MT_daughters != nullptr ) return;
315
316 G4int index;
317
318#ifdef G4VERBOSE
319 if (verboseLevel>1) G4cout << "G4VDecayChannel::FillDaughters()" << G4endl;
320#endif
321 if (G4MT_daughters != nullptr)
322 {
323 delete [] G4MT_daughters;
324 G4MT_daughters = nullptr;
325 }
326
327 // parent mass
329 G4double parentmass = G4MT_parent->GetPDGMass();
330
331 //
332 G4double sumofdaughtermass = 0.0;
333 G4double sumofdaughterwidthsq = 0.0;
334
335 if ((numberOfDaughters <=0) || (daughters_name == nullptr) )
336 {
337#ifdef G4VERBOSE
338 if (verboseLevel>0)
339 {
340 G4cout << "G4VDecayChannel::FillDaughters() - "
341 << "[ " << G4MT_parent->GetParticleName() << " ]"
342 << "numberOfDaughters is not defined yet";
343 }
344#endif
345 G4MT_daughters = nullptr;
346 G4Exception("G4VDecayChannel::FillDaughters()",
347 "PART011", FatalException,
348 "Cannot fill daughters: numberOfDaughters is not defined yet");
349 }
350
351 // create and set the array of pointers to daughter particles
353 if (G4MT_daughters_mass != nullptr) delete [] G4MT_daughters_mass;
354 if (G4MT_daughters_width != nullptr) delete [] G4MT_daughters_width;
357 // loop over all daughters
358 for (index=0; index<numberOfDaughters; ++index)
359 {
360 if (daughters_name[index] == nullptr)
361 {
362 // daughter name is not defined
363#ifdef G4VERBOSE
364 if (verboseLevel>0)
365 {
366 G4cout << "G4VDecayChannel::FillDaughters() - "
367 << "[ " << G4MT_parent->GetParticleName() << " ]"
368 << index << "-th daughter is not defined yet" << G4endl;
369 }
370#endif
371 G4MT_daughters[index] = nullptr;
372 G4Exception("G4VDecayChannel::FillDaughters()",
373 "PART011", FatalException,
374 "Cannot fill daughters: name of daughter is not defined yet");
375 }
376 // search daughter particles in the particle table
378 if (G4MT_daughters[index] == nullptr )
379 {
380 // cannot find the daughter particle
381#ifdef G4VERBOSE
382 if (verboseLevel>0)
383 {
384 G4cout << "G4VDecayChannel::FillDaughters() - "
385 << "[ " << G4MT_parent->GetParticleName() << " ]"
386 << index << ":" << *daughters_name[index]
387 << " is not defined !!" << G4endl;
388 G4cout << " The BR of this decay mode is set to zero." << G4endl;
389 }
390#endif
391 SetBR(0.0);
392 return;
393 }
394#ifdef G4VERBOSE
395 if (verboseLevel>1)
396 {
397 G4cout << index << ":" << *daughters_name[index];
398 G4cout << ":" << G4MT_daughters[index] << G4endl;
399 }
400#endif
402 G4double d_width = G4MT_daughters[index]->GetPDGWidth();
403 G4MT_daughters_width[index] = d_width;
404 sumofdaughtermass += G4MT_daughters[index]->GetPDGMass();
405 sumofdaughterwidthsq += d_width*d_width;
406 } // end loop over all daughters
407
408 // check sum of daghter mass
409 G4double widthMass = std::sqrt(G4MT_parent->GetPDGWidth()
411 +sumofdaughterwidthsq);
412 if ( (G4MT_parent->GetParticleType() != "nucleus")
413 && (numberOfDaughters !=1)
414 && (sumofdaughtermass > parentmass + rangeMass*widthMass) )
415 {
416 // !!! illegal mass !!!
417#ifdef G4VERBOSE
418 if (GetVerboseLevel()>0)
419 {
420 G4cout << "G4VDecayChannel::FillDaughters() - "
421 << "[ " << G4MT_parent->GetParticleName() << " ]"
422 << " Energy/Momentum conserevation breaks " << G4endl;
423 if (GetVerboseLevel()>1)
424 {
425 G4cout << " parent:" << *parent_name
426 << " mass:" << parentmass/GeV << "[GeV/c/c]" << G4endl;
427 for (index=0; index<numberOfDaughters; ++index)
428 {
429 G4cout << " daughter " << index << ":" << *daughters_name[index]
430 << " mass:" << G4MT_daughters[index]->GetPDGMass()/GeV
431 << "[GeV/c/c]" << G4endl;
432 }
433 }
434 }
435#endif
436 }
437}
@ FatalException
static constexpr double GeV
Definition: G4SIunits.hh:203
const G4String & GetParticleType() const
G4double GetPDGWidth() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void SetBR(G4double value)
G4ParticleTable * particletable

References G4VDecayChannel::CheckAndFillParent(), G4VDecayChannel::daughters_name, G4VDecayChannel::daughtersMutex, FatalException, G4ParticleTable::FindParticle(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::G4MT_parent, G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGMass(), G4ParticleDefinition::GetPDGWidth(), G4VDecayChannel::GetVerboseLevel(), GeV, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, G4VDecayChannel::particletable, G4VDecayChannel::rangeMass, G4VDecayChannel::SetBR(), and G4VDecayChannel::verboseLevel.

Referenced by G4VDecayChannel::CheckAndFillDaughters().

◆ FillParent()

void G4VDecayChannel::FillParent ( )
privateinherited

Definition at line 440 of file G4VDecayChannel.cc.

441{
442 G4AutoLock lock(&parentMutex);
443
444 // Double check, check again if another thread has already filled this, in
445 // case do not need to do anything
446 if ( G4MT_parent != nullptr ) return;
447
448 if (parent_name == nullptr)
449 {
450 // parent name is not defined
451#ifdef G4VERBOSE
452 if (verboseLevel>0)
453 {
454 G4cout << "G4VDecayChannel::FillParent() - "
455 << "parent name is not defined !!" << G4endl;
456 }
457#endif
458 G4MT_parent = nullptr;
459 G4Exception("G4VDecayChannel::FillParent()",
460 "PART012", FatalException,
461 "Cannot fill parent: parent name is not defined yet");
462 return;
463 }
464
465 // search parent particle in the particle table
467 if (G4MT_parent == nullptr)
468 {
469 // parent particle does not exist
470#ifdef G4VERBOSE
471 if (verboseLevel>0)
472 {
473 G4cout << "G4VDecayChannel::FillParent() - "
474 << *parent_name << " does not exist !!" << G4endl;
475 }
476#endif
477 G4Exception("G4VDecayChannel::FillParent()",
478 "PART012", FatalException,
479 "Cannot fill parent: parent does not exist");
480 return;
481 }
483}
G4double G4MT_parent_mass

References FatalException, G4ParticleTable::FindParticle(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_parent, G4VDecayChannel::G4MT_parent_mass, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::parent_name, G4VDecayChannel::parentMutex, G4VDecayChannel::particletable, and G4VDecayChannel::verboseLevel.

Referenced by G4VDecayChannel::CheckAndFillParent().

◆ GetAngularMomentum()

G4int G4VDecayChannel::GetAngularMomentum ( )
inherited

Definition at line 492 of file G4VDecayChannel.cc.

493{
494 // determine angular momentum
495
496 // fill pointers to daughter particles if not yet set
498
499 const G4int PiSpin = G4MT_parent->GetPDGiSpin();
500 const G4int PParity = G4MT_parent->GetPDGiParity();
501 if (2==numberOfDaughters) // up to now we can only handle two particle decays
502 {
503 const G4int D1iSpin = G4MT_daughters[0]->GetPDGiSpin();
504 const G4int D1Parity = G4MT_daughters[0]->GetPDGiParity();
505 const G4int D2iSpin = G4MT_daughters[1]->GetPDGiSpin();
506 const G4int D2Parity = G4MT_daughters[1]->GetPDGiParity();
507 const G4int MiniSpin = std::abs (D1iSpin - D2iSpin);
508 const G4int MaxiSpin = D1iSpin + D2iSpin;
509 const G4int lMax = (PiSpin+D1iSpin+D2iSpin)/2; // l is always int
510 G4int lMin;
511#ifdef G4VERBOSE
512 if (verboseLevel>1)
513 {
514 G4cout << "iSpin: " << PiSpin << " -> " << D1iSpin
515 << " + " << D2iSpin << G4endl;
516 G4cout << "2*jmin, 2*jmax, lmax " << MiniSpin << " " << MaxiSpin
517 << " " << lMax << G4endl;
518 }
519#endif
520 for (G4int j=MiniSpin; j<=MaxiSpin; j+=2) // loop over all possible
521 { // spin couplings
522 lMin = std::abs(PiSpin-j)/2;
523#ifdef G4VERBOSE
524 if (verboseLevel>1)
525 G4cout << "-> checking 2*j=" << j << G4endl;
526#endif
527 for (G4int l=lMin; l<=lMax; ++l)
528 {
529#ifdef G4VERBOSE
530 if (verboseLevel>1)
531 G4cout << " checking l=" << l << G4endl;
532#endif
533 if (l%2==0)
534 {
535 if (PParity == D1Parity*D2Parity) // check parity for this l
536 return l;
537 }
538 else
539 {
540 if (PParity == -1*D1Parity*D2Parity) // check parity for this l
541 return l;
542 }
543 }
544 }
545 }
546 else
547 {
548 G4Exception("G4VDecayChannel::GetAngularMomentum()",
549 "PART111", JustWarning,
550 "Sorry, can't handle 3 particle decays (up to now)");
551 return 0;
552 }
553 G4Exception ("G4VDecayChannel::GetAngularMomentum()",
554 "PART111", JustWarning,
555 "Can't find angular momentum for this decay");
556 return 0;
557}

References G4VDecayChannel::CheckAndFillDaughters(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4ParticleDefinition::GetPDGiParity(), G4ParticleDefinition::GetPDGiSpin(), JustWarning, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

◆ GetBR()

G4double G4VDecayChannel::GetBR ( ) const
inlineinherited

◆ GetDaughter()

G4ParticleDefinition * G4VDecayChannel::GetDaughter ( G4int  anIndex)
inlineinherited

Definition at line 201 of file G4VDecayChannel.hh.

202{
203 // pointers to daughter particles are filled, if they are not set yet
205
206 // get the pointer to a daughter particle
207 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
208 {
209 return G4MT_daughters[anIndex];
210 }
211 else
212 {
213 if (verboseLevel>0)
214 G4cout << "G4VDecayChannel::GetDaughter index out of range "
215 << anIndex << G4endl;
216 return nullptr;
217 }
218}

References G4VDecayChannel::CheckAndFillDaughters(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

Referenced by G4IonTable::CreateIon(), G4KineticTrack::Decay(), G4KineticTrack::G4KineticTrack(), G4HtmlPPReporter::GeneratePropertyTable(), G4TextPPReporter::GeneratePropertyTable(), G4NuclearDecay::GetDaughterNucleus(), G4MTRunManagerKernel::SetUpDecayChannels(), and G4TaskRunManagerKernel::SetUpDecayChannels().

◆ GetDaughterMass()

G4double G4VDecayChannel::GetDaughterMass ( G4int  anIndex) const
inlineinherited

Definition at line 239 of file G4VDecayChannel.hh.

240{
241 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
242 {
243 return G4MT_daughters_mass[anIndex];
244 }
245 else
246 {
247 if (verboseLevel>0)
248 {
249 G4cout << "G4VDecayChannel::GetDaughterMass ";
250 G4cout << "index out of range " << anIndex << G4endl;
251 }
252 return 0.0;
253 }
254}

References G4cout, G4endl, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

◆ GetDaughterName()

const G4String & G4VDecayChannel::GetDaughterName ( G4int  anIndex) const
inlineinherited

◆ GetKinematicsName()

const G4String & G4VDecayChannel::GetKinematicsName ( ) const
inlineinherited

◆ GetNoName()

const G4String & G4VDecayChannel::GetNoName ( ) const
privateinherited

Definition at line 579 of file G4VDecayChannel.cc.

580{
581 return noName;
582}
static const G4String noName

References G4VDecayChannel::noName.

Referenced by G4VDecayChannel::GetDaughterName().

◆ GetNumberOfDaughters()

G4int G4VDecayChannel::GetNumberOfDaughters ( ) const
inlineinherited

◆ GetParent()

G4ParticleDefinition * G4VDecayChannel::GetParent ( )
inlineinherited

Definition at line 257 of file G4VDecayChannel.hh.

258{
259 // the pointer to the parent particle is filled, if it is not set yet
261 // get the pointer to the parent particle
262 return G4MT_parent;
263}

References G4VDecayChannel::CheckAndFillParent(), and G4VDecayChannel::G4MT_parent.

Referenced by G4DecayTable::Insert().

◆ GetParentMass()

G4double G4VDecayChannel::GetParentMass ( ) const
inlineinherited

Definition at line 272 of file G4VDecayChannel.hh.

273{
274 return G4MT_parent_mass;
275}

References G4VDecayChannel::G4MT_parent_mass.

◆ GetParentName()

const G4String & G4VDecayChannel::GetParentName ( ) const
inlineinherited

◆ GetPolarization()

const G4ThreeVector & G4VDecayChannel::GetPolarization ( ) const
inlineinherited

Definition at line 331 of file G4VDecayChannel.hh.

332{
333 return parent_polarization;
334}

References G4VDecayChannel::parent_polarization.

◆ GetRangeMass()

G4double G4VDecayChannel::GetRangeMass ( ) const
inlineinherited

Definition at line 316 of file G4VDecayChannel.hh.

317{
318 return rangeMass;
319}

References G4VDecayChannel::rangeMass.

◆ GetVerboseLevel()

G4int G4VDecayChannel::GetVerboseLevel ( ) const
inlineinherited

◆ IsOKWithParentMass()

G4bool G4VDecayChannel::IsOKWithParentMass ( G4double  parentMass)
virtualinherited

Reimplemented in G4PhaseSpaceDecayChannel.

Definition at line 607 of file G4VDecayChannel.cc.

608{
609 G4double sumOfDaughterMassMin = 0.0;
612 // skip one body decay
613 if (numberOfDaughters==1) return true;
614
615 for (G4int index=0; index<numberOfDaughters; ++index)
616 {
617 sumOfDaughterMassMin +=
619 }
620 return (parentMass >= sumOfDaughterMassMin);
621}

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::rangeMass.

Referenced by G4Decay::DecayIt(), G4MuonicAtomDecay::DecayIt(), and G4PhaseSpaceDecayChannel::IsOKWithParentMass().

◆ operator!=()

G4bool G4VDecayChannel::operator!= ( const G4VDecayChannel r) const
inlineinherited

Definition at line 69 of file G4VDecayChannel.hh.

69{ return (this != &r); }

◆ operator<()

G4bool G4VDecayChannel::operator< ( const G4VDecayChannel right) const
inlineinherited

Definition at line 195 of file G4VDecayChannel.hh.

196{
197 return (this->rbranch < right.rbranch);
198}

References G4VDecayChannel::rbranch.

◆ operator=()

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

Definition at line 68 of file G4MuonDecayChannelWithSpin.cc.

70{
71 if (this != &right)
72 {
75 rbranch = right.rbranch;
76
77 // copy parent name
78 delete parent_name;
79 parent_name = new G4String(*right.parent_name);
80
81 // clear daughters_name array
83
84 // recreate array
86 if ( numberOfDaughters > 0 )
87 {
89 // copy daughters name
90 for (G4int index=0; index<numberOfDaughters; ++index)
91 {
92 daughters_name[index] = new G4String(*right.daughters_name[index]);
93 }
94 }
95 }
96 return *this;
97}

References G4VDecayChannel::ClearDaughtersName(), G4VDecayChannel::daughters_name, G4VDecayChannel::kinematics_name, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, G4VDecayChannel::rbranch, and G4VDecayChannel::verboseLevel.

◆ operator==()

G4bool G4VDecayChannel::operator== ( const G4VDecayChannel r) const
inlineinherited

Definition at line 68 of file G4VDecayChannel.hh.

68{ return (this == &r); }

◆ R_c()

G4double G4MuonDecayChannelWithSpin::R_c ( G4double  x,
G4double  omega 
)
private

Definition at line 304 of file G4MuonDecayChannelWithSpin.cc.

305{
306 G4int n_max = (G4int)(100.*x);
307
308 if(n_max<10)n_max=10;
309
310 G4double L2 = 0.0;
311
312 for(G4int n=1; n<=n_max; ++n)
313 {
314 L2 += std::pow(x,n)/(n*n);
315 }
316
317 G4double r_c;
318
319 r_c = 2.*L2-(pi*pi/3.)-2.;
320 r_c = r_c + omega * (1.5+2.*std::log((1.-x)/x));
321 r_c = r_c - std::log(x)*(2.*std::log(x)-1.);
322 r_c = r_c + (3.*std::log(x)-1.-1./x)*std::log(1.-x);
323
324 return r_c;
325}
static constexpr double pi
Definition: G4SIunits.hh:55

References CLHEP::detail::n, and pi.

Referenced by F_c(), and F_theta().

◆ SetBR()

void G4VDecayChannel::SetBR ( G4double  value)
inherited

◆ SetDaughter() [1/2]

void G4VDecayChannel::SetDaughter ( G4int  anIndex,
const G4ParticleDefinition particle_type 
)
inherited

◆ SetDaughter() [2/2]

void G4VDecayChannel::SetDaughter ( G4int  anIndex,
const G4String particle_name 
)
inherited

Definition at line 234 of file G4VDecayChannel.cc.

235{
236 // check numberOfDaughters is positive
237 if (numberOfDaughters<=0)
238 {
239#ifdef G4VERBOSE
240 if (verboseLevel>0)
241 {
242 G4cout << "G4VDecayChannel::SetDaughter() - "
243 << "Number of daughters is not defined" << G4endl;
244 }
245#endif
246 return;
247 }
248
249 // An analysis of this code, shows that this method is called
250 // only in the constructor of derived classes.
251 // The general idea of this method is probably to support
252 // the possibility to re-define daughters on the fly, however
253 // this design is extremely problematic for MT mode, we thus
254 // require (as practically happens) that the method is called only
255 // at construction, i.e. when G4MT_daughters == 0
256 // moreover this method can be called only after SetNumberOfDaugthers()
257 // has been called (see previous if), in such a case daughters_name != nullptr
258 //
259 if ( daughters_name == nullptr )
260 {
261 G4Exception("G4VDecayChannel::SetDaughter()", "PART112", FatalException,
262 "Trying to add a daughter without specifying number of secondaries!");
263 return;
264 }
265 if ( G4MT_daughters != nullptr )
266 {
267 G4Exception("G4VDecayChannel::SetDaughter()", "PART111", FatalException,
268 "Trying to modify a daughter of a decay channel, \
269 but decay channel already has daughters.");
270 return;
271 }
272
273 // check an index
274 if ( (anIndex<0) || (anIndex>=numberOfDaughters) )
275 {
276#ifdef G4VERBOSE
277 if (verboseLevel>0)
278 {
279 G4cout << "G4VDecayChannel::SetDaughter() - "
280 << "index out of range " << anIndex << G4endl;
281 }
282#endif
283 }
284 else
285 {
286 // fill the name
287 daughters_name[anIndex] = new G4String(particle_name);
288#ifdef G4VERBOSE
289 if (verboseLevel>1)
290 {
291 G4cout << "G4VDecayChannel::SetDaughter[" << anIndex <<"] :";
292 G4cout << daughters_name[anIndex] << ":"
293 << *daughters_name[anIndex]<<G4endl;
294 }
295#endif
296 }
297}

References G4VDecayChannel::daughters_name, FatalException, G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

◆ SetNumberOfDaughters()

void G4VDecayChannel::SetNumberOfDaughters ( G4int  value)
inherited

◆ SetParent() [1/2]

void G4VDecayChannel::SetParent ( const G4ParticleDefinition particle_type)
inherited

◆ SetParent() [2/2]

void G4VDecayChannel::SetParent ( const G4String particle_name)
inlineinherited

Definition at line 278 of file G4VDecayChannel.hh.

279{
280 if (parent_name != nullptr) delete parent_name;
281 parent_name = new G4String(particle_name);
282 G4MT_parent = nullptr;
283}

References G4VDecayChannel::G4MT_parent, and G4VDecayChannel::parent_name.

◆ SetPolarization()

void G4VDecayChannel::SetPolarization ( const G4ThreeVector polar)
inlineinherited

◆ SetRangeMass()

void G4VDecayChannel::SetRangeMass ( G4double  val)
inlineinherited

Definition at line 322 of file G4VDecayChannel.hh.

322{ if(val>=0.) rangeMass=val; }

References G4VDecayChannel::rangeMass.

◆ SetVerboseLevel()

void G4VDecayChannel::SetVerboseLevel ( G4int  value)
inlineinherited

Definition at line 304 of file G4VDecayChannel.hh.

305{
306 verboseLevel = value;
307}

References G4VDecayChannel::verboseLevel.

Referenced by G4Decay::DecayIt(), and G4MuonicAtomDecay::DecayIt().

Field Documentation

◆ daughters_name

G4String** G4VDecayChannel::daughters_name = nullptr
protectedinherited

◆ daughtersMutex

G4Mutex G4VDecayChannel::daughtersMutex
protectedinherited

◆ G4MT_daughters

G4ParticleDefinition** G4VDecayChannel::G4MT_daughters = nullptr
protectedinherited

◆ G4MT_daughters_mass

G4double* G4VDecayChannel::G4MT_daughters_mass = nullptr
protectedinherited

◆ G4MT_daughters_width

G4double* G4VDecayChannel::G4MT_daughters_width = nullptr
protectedinherited

◆ G4MT_parent

G4ParticleDefinition* G4VDecayChannel::G4MT_parent = nullptr
protectedinherited

◆ G4MT_parent_mass

G4double G4VDecayChannel::G4MT_parent_mass = 0.0
protectedinherited

◆ kinematics_name

G4String G4VDecayChannel::kinematics_name = ""
protectedinherited

◆ noName

const G4String G4VDecayChannel::noName = " "
staticprotectedinherited

Definition at line 170 of file G4VDecayChannel.hh.

Referenced by G4VDecayChannel::GetNoName().

◆ numberOfDaughters

G4int G4VDecayChannel::numberOfDaughters = 0
protectedinherited

◆ parent_name

G4String* G4VDecayChannel::parent_name = nullptr
protectedinherited

◆ parent_polarization

G4ThreeVector G4VDecayChannel::parent_polarization
protectedinherited

◆ parentMutex

G4Mutex G4VDecayChannel::parentMutex
protectedinherited

◆ particletable

G4ParticleTable* G4VDecayChannel::particletable = nullptr
protectedinherited

◆ rangeMass

G4double G4VDecayChannel::rangeMass = 2.5
protectedinherited

◆ rbranch

G4double G4VDecayChannel::rbranch = 0.0
protectedinherited

◆ verboseLevel

G4int G4VDecayChannel::verboseLevel = 1
protectedinherited

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