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

#include <G4BetaPlusDecay.hh>

Inheritance diagram for G4BetaPlusDecay:
G4NuclearDecay G4VDecayChannel

Public Member Functions

virtual G4DecayProductsDecayIt (G4double)
 
void DumpInfo ()
 
virtual void DumpNuclearInfo ()
 
 G4BetaPlusDecay (const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &endpointE, const G4double &ex, const G4Ions::G4FloatLevelBase &flb, const G4BetaDecayType &type)
 
G4int GetAngularMomentum ()
 
G4double GetBR () const
 
G4ParticleDefinitionGetDaughter (G4int anIndex)
 
G4double GetDaughterExcitation ()
 
G4double GetDaughterMass (G4int anIndex) const
 
const G4StringGetDaughterName (G4int anIndex) const
 
G4ParticleDefinitionGetDaughterNucleus ()
 
G4RadioactiveDecayMode GetDecayMode ()
 
G4Ions::G4FloatLevelBase GetFloatingLevel ()
 
G4double GetHLThreshold ()
 
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 SetHLThreshold (G4double HLT)
 
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 ~G4BetaPlusDecay ()
 

Protected Member Functions

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

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
 
const G4RadioactiveDecayMode theMode
 
G4int verboseLevel = 1
 

Static Protected Attributes

static const G4String noName = " "
 

Private Member Functions

void FillDaughters ()
 
void FillParent ()
 
const G4StringGetNoName () const
 
void SetUpBetaSpectrumSampler (const G4int &parentZ, const G4int &parentA, const G4BetaDecayType &type)
 

Private Attributes

const G4double daughterEx
 
const G4double endpointEnergy
 
const G4Ions::G4FloatLevelBase floatingLevel
 
G4double halflifeThreshold
 
G4RandGeneralspectrumSampler
 

Detailed Description

Definition at line 44 of file G4BetaPlusDecay.hh.

Constructor & Destructor Documentation

◆ G4BetaPlusDecay()

G4BetaPlusDecay::G4BetaPlusDecay ( const G4ParticleDefinition theParentNucleus,
const G4double theBR,
const G4double endpointE,
const G4double ex,
const G4Ions::G4FloatLevelBase flb,
const G4BetaDecayType type 
)

Definition at line 45 of file G4BetaPlusDecay.cc.

50 : G4NuclearDecay("beta+ decay", BetaPlus, excitationE, flb),
52{
53 SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
54 SetBR(branch);
55
57 G4IonTable* theIonTable =
59 G4int daughterZ = theParentNucleus->GetAtomicNumber() - 1;
60 G4int daughterA = theParentNucleus->GetAtomicMass();
61 SetDaughter(0, theIonTable->GetIon(daughterZ, daughterA, excitationE, flb) );
62 SetUpBetaSpectrumSampler(daughterZ, daughterA, betaType);
63 SetDaughter(1, "e+");
64 SetDaughter(2, "nu_e");
65}
int G4int
Definition: G4Types.hh:85
void SetUpBetaSpectrumSampler(const G4int &parentZ, const G4int &parentA, const G4BetaDecayType &type)
const G4double endpointEnergy
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
G4NuclearDecay(const G4String &channelName, const G4RadioactiveDecayMode &mode, const G4double &excitation, const G4Ions::G4FloatLevelBase &floatingLevel)
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
void SetBR(G4double value)
void SetNumberOfDaughters(G4int value)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
void SetParent(const G4ParticleDefinition *particle_type)
static constexpr double electron_mass_c2

References G4ParticleDefinition::GetAtomicMass(), G4ParticleDefinition::GetAtomicNumber(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4ParticleTable::GetParticleTable(), G4VDecayChannel::SetBR(), G4VDecayChannel::SetDaughter(), G4VDecayChannel::SetNumberOfDaughters(), G4VDecayChannel::SetParent(), and SetUpBetaSpectrumSampler().

◆ ~G4BetaPlusDecay()

G4BetaPlusDecay::~G4BetaPlusDecay ( )
virtual

Definition at line 68 of file G4BetaPlusDecay.cc.

69{
70 delete spectrumSampler;
71}
G4RandGeneral * spectrumSampler

References spectrumSampler.

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}
#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=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4NeutronBetaDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4TauLeptonicDecayChannel::operator=(), G4VDecayChannel::operator=(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::~G4VDecayChannel().

◆ DecayIt()

G4DecayProducts * G4BetaPlusDecay::DecayIt ( G4double  )
virtual

Implements G4VDecayChannel.

Definition at line 74 of file G4BetaPlusDecay.cc.

75{
76 // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)
78
79 // Fill G4MT_daughters with e-, nu and residual nucleus (stored by SetDaughter)
81
82 G4double parentMass = G4MT_parent->GetPDGMass();
84 G4double nucleusMass = G4MT_daughters[0]->GetPDGMass();
85 // Set up final state
86 // parentParticle is set at rest here because boost with correct momentum
87 // is done later
88 G4DynamicParticle parentParticle(G4MT_parent, G4ThreeVector(0,0,0), 0.0);
89 G4DecayProducts* products = new G4DecayProducts(parentParticle);
90
91 if (spectrumSampler) {
92 // Generate positron isotropic in angle, with energy from stored spectrum
93 G4double eKE = endpointEnergy*spectrumSampler->shoot(G4Random::getTheEngine() );
94 G4double eMomentum = std::sqrt(eKE*(eKE + 2.*eMass) );
95
96 G4double cosTheta = 2.*G4UniformRand() - 1.0;
97 G4double sinTheta = std::sqrt(1.0 - cosTheta*cosTheta);
99 G4double sinPhi = std::sin(phi);
100 G4double cosPhi = std::cos(phi);
101
102 G4ParticleMomentum eDirection(sinTheta*cosPhi, sinTheta*sinPhi, cosTheta);
103 G4DynamicParticle* dynamicPositron
104 = new G4DynamicParticle(G4MT_daughters[1], eDirection*eMomentum);
105 products->PushProducts(dynamicPositron);
106
107 // Generate neutrino with angle relative to positron, and energy from
108 // energy-momentum conservation using endpoint energy of reaction
109 G4double cosThetaENu = 2.*G4UniformRand() - 1.;
110 G4double eTE = eMass + eKE;
111 G4double nuEnergy = ((endpointEnergy - eKE)*(parentMass + nucleusMass - eTE)
112 - eMomentum*eMomentum)/(parentMass - eTE + eMomentum*cosThetaENu)/2.;
113
114 G4double sinThetaENu = std::sqrt(1.0 - cosThetaENu*cosThetaENu);
115 phi = twopi*G4UniformRand()*rad;
116 G4double sinPhiNu = std::sin(phi);
117 G4double cosPhiNu = std::cos(phi);
118
119 G4ParticleMomentum nuDirection;
120 nuDirection.setX(sinThetaENu*cosPhiNu*cosTheta*cosPhi -
121 sinThetaENu*sinPhiNu*sinPhi + cosThetaENu*sinTheta*cosPhi);
122 nuDirection.setY(sinThetaENu*cosPhiNu*cosTheta*sinPhi +
123 sinThetaENu*sinPhiNu*cosPhi + cosThetaENu*sinTheta*sinPhi);
124 nuDirection.setZ(-sinThetaENu*cosPhiNu*sinTheta + cosThetaENu*cosTheta);
125
126 G4DynamicParticle* dynamicNeutrino
127 = new G4DynamicParticle(G4MT_daughters[2], nuDirection*nuEnergy);
128 products->PushProducts(dynamicNeutrino);
129
130 // Generate daughter nucleus from sum of positron and neutrino 4-vectors:
131 // p_D = - p_e - p_nu
132 G4DynamicParticle* dynamicDaughter =
134 -eDirection*eMomentum - nuDirection*nuEnergy);
135 products->PushProducts(dynamicDaughter);
136
137 } else {
138 // positron energy below threshold -> no decay
139 G4DynamicParticle* noDecay =
141 products->PushProducts(noDecay);
142 }
143
144 // Check energy conservation against endpoint value, not nuclear masses
145 /*
146 G4int nProd = products->entries();
147 G4DynamicParticle* temp = 0;
148 G4double Esum = 0.0;
149 for (G4int i = 0; i < nProd; i++) {
150 temp = products->operator[](i);
151 Esum += temp->GetKineticEnergy();
152 }
153 G4double eCons = (endpointEnergy - Esum)/keV;
154 if (eCons > 0.001) G4cout << " Beta+ check: eCons (keV) = " << eCons << G4endl;
155 */
156 return products;
157}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double rad
Definition: G4SIunits.hh:129
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52
void setY(double)
void setZ(double)
void setX(double)
G4int PushProducts(G4DynamicParticle *aParticle)
void CheckAndFillDaughters()

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), endpointEnergy, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4UniformRand, G4ParticleDefinition::GetPDGMass(), G4DecayProducts::PushProducts(), rad, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), spectrumSampler, and twopi.

◆ 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().

◆ DumpNuclearInfo()

void G4BetaPlusDecay::DumpNuclearInfo ( )
virtual

Implements G4NuclearDecay.

Definition at line 197 of file G4BetaPlusDecay.cc.

198{
199 G4cout << " G4BetaPlusDecay for parent nucleus " << GetParentName() << G4endl;
200 G4cout << " decays to " << GetDaughterName(0) << " , " << GetDaughterName(1)
201 << " and " << GetDaughterName(2) << " with branching ratio " << GetBR()
202 << "% and endpoint energy " << endpointEnergy/keV << " keV " << G4endl;
203}
static constexpr double keV
Definition: G4SIunits.hh:202
G4double GetBR() const
const G4String & GetParentName() const
const G4String & GetDaughterName(G4int anIndex) const

References endpointEnergy, G4cout, G4endl, G4VDecayChannel::GetBR(), G4VDecayChannel::GetDaughterName(), G4VDecayChannel::GetParentName(), and keV.

◆ 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().

◆ 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
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double GeV
Definition: G4SIunits.hh:203
const G4String & GetParticleType() const
G4double GetPDGWidth() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4int GetVerboseLevel() const
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}
@ JustWarning

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().

◆ GetDaughterExcitation()

G4double G4NuclearDecay::GetDaughterExcitation ( )
inlineinherited

Definition at line 55 of file G4NuclearDecay.hh.

55{return daughterEx;}
const G4double daughterEx

References G4NuclearDecay::daughterEx.

Referenced by G4Radioactivation::CalculateChainsFromParent().

◆ 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

◆ GetDaughterNucleus()

G4ParticleDefinition * G4NuclearDecay::GetDaughterNucleus ( )
inlineinherited

Definition at line 59 of file G4NuclearDecay.hh.

59{return GetDaughter(0);}
G4ParticleDefinition * GetDaughter(G4int anIndex)

References G4VDecayChannel::GetDaughter().

Referenced by G4Radioactivation::CalculateChainsFromParent().

◆ GetDecayMode()

G4RadioactiveDecayMode G4NuclearDecay::GetDecayMode ( )
inlineinherited

Definition at line 53 of file G4NuclearDecay.hh.

53{return theMode;}
const G4RadioactiveDecayMode theMode

References G4NuclearDecay::theMode.

Referenced by G4Radioactivation::CalculateChainsFromParent(), and G4RadioactiveDecay::LoadDecayTable().

◆ GetFloatingLevel()

G4Ions::G4FloatLevelBase G4NuclearDecay::GetFloatingLevel ( )
inlineinherited

Definition at line 57 of file G4NuclearDecay.hh.

57{return floatingLevel;}
const G4Ions::G4FloatLevelBase floatingLevel

References G4NuclearDecay::floatingLevel.

◆ GetHLThreshold()

G4double G4NuclearDecay::GetHLThreshold ( )
inlineinherited

Definition at line 62 of file G4NuclearDecay.hh.

62{return halflifeThreshold;}
G4double halflifeThreshold

References G4NuclearDecay::halflifeThreshold.

◆ 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}
G4ThreeVector parent_polarization

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==()

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

Definition at line 68 of file G4VDecayChannel.hh.

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

◆ 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.

◆ SetHLThreshold()

void G4NuclearDecay::SetHLThreshold ( G4double  HLT)
inlineinherited

Definition at line 61 of file G4NuclearDecay.hh.

61{halflifeThreshold = HLT;}

References G4NuclearDecay::halflifeThreshold.

◆ 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.

◆ SetUpBetaSpectrumSampler()

void G4BetaPlusDecay::SetUpBetaSpectrumSampler ( const G4int parentZ,
const G4int parentA,
const G4BetaDecayType type 
)
private

Definition at line 161 of file G4BetaPlusDecay.cc.

164{
166 G4BetaDecayCorrections corrections(-daughterZ, daughterA);
167 spectrumSampler = 0;
168
169 // Check for cases in which Q < 2Me (e.g. z67.a162)
170 if (e0 > 0.) {
171 // Array to store spectrum pdf
172 G4int npti = 100;
173 G4double* pdf = new G4double[npti];
174
175 G4double e; // Total positron energy in units of electron mass
176 G4double p; // Positron momentum in units of electron mass
177 G4double f; // Spectral shap function
178 for (G4int ptn = 0; ptn < npti; ptn++) {
179 // Calculate simple phase space
180 e = 1. + e0*(ptn + 0.5)/G4double(npti);
181 p = std::sqrt(e*e - 1.);
182 f = p*e*(e0 - e + 1.)*(e0 - e + 1.);
183
184 // Apply Fermi factor to get allowed shape
185 f *= corrections.FermiFunction(e);
186
187 // Apply shape factor for forbidden transitions
188 f *= corrections.ShapeFactor(betaType, p, e0-e+1.);
189 pdf[ptn] = f;
190 }
191 spectrumSampler = new G4RandGeneral(pdf, npti);
192 delete[] pdf;
193 }
194}
#define G4RandGeneral
Definition: Randomize.hh:49

References CLHEP::electron_mass_c2, endpointEnergy, G4BetaDecayCorrections::FermiFunction(), G4RandGeneral, G4BetaDecayCorrections::ShapeFactor(), and spectrumSampler.

Referenced by G4BetaPlusDecay().

◆ 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

◆ daughterEx

const G4double G4NuclearDecay::daughterEx
privateinherited

Definition at line 71 of file G4NuclearDecay.hh.

Referenced by G4NuclearDecay::GetDaughterExcitation().

◆ daughters_name

G4String** G4VDecayChannel::daughters_name = nullptr
protectedinherited

◆ daughtersMutex

G4Mutex G4VDecayChannel::daughtersMutex
protectedinherited

◆ endpointEnergy

const G4double G4BetaPlusDecay::endpointEnergy
private

Definition at line 62 of file G4BetaPlusDecay.hh.

Referenced by DecayIt(), DumpNuclearInfo(), and SetUpBetaSpectrumSampler().

◆ floatingLevel

const G4Ions::G4FloatLevelBase G4NuclearDecay::floatingLevel
privateinherited

Definition at line 72 of file G4NuclearDecay.hh.

Referenced by G4NuclearDecay::GetFloatingLevel().

◆ 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

◆ halflifeThreshold

G4double G4NuclearDecay::halflifeThreshold
privateinherited

◆ 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

◆ spectrumSampler

G4RandGeneral* G4BetaPlusDecay::spectrumSampler
private

Definition at line 63 of file G4BetaPlusDecay.hh.

Referenced by DecayIt(), SetUpBetaSpectrumSampler(), and ~G4BetaPlusDecay().

◆ theMode

const G4RadioactiveDecayMode G4NuclearDecay::theMode
protectedinherited

◆ verboseLevel

G4int G4VDecayChannel::verboseLevel = 1
protectedinherited

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