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

#include <G4MuonicAtomDecay.hh>

+ Inheritance diagram for G4MuonicAtomDecay:

Public Member Functions

 G4MuonicAtomDecay (G4HadronicInteraction *hiptr=nullptr, const G4String &processName="MuonicAtomDecay")
 
virtual ~G4MuonicAtomDecay ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &theTrack, const G4Step &theStep)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &theTrack, const G4Step &theStep)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 50 of file G4MuonicAtomDecay.hh.

Constructor & Destructor Documentation

◆ G4MuonicAtomDecay()

G4MuonicAtomDecay::G4MuonicAtomDecay ( G4HadronicInteraction hiptr = nullptr,
const G4String processName = "MuonicAtomDecay" 
)
explicit

Definition at line 67 of file G4MuonicAtomDecay.cc.

70 fMuMass(G4MuonMinus::MuonMinus()->GetPDGMass()),
71 cmptr(hiptr),
72 verboseLevel(0)
73{
74 // This is not a hadronic process; assume it is a kind of decay
75 enableAtRestDoIt = true;
76 enablePostStepDoIt = true; // it is a streach; fixme
77 theProcessSubType = 221; // (see G4DecayProcessType.hh) fixme
78 if (!cmptr) {
79 // cmptr = new G4CascadeInterface(); // Bertini - Pointer owned by InteractionRegistry
80 cmptr = new G4MuMinusCapturePrecompound(); // Precompound
81 }
82}
@ fDecay
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:365
G4int theProcessSubType
Definition: G4VProcess.hh:353

◆ ~G4MuonicAtomDecay()

G4MuonicAtomDecay::~G4MuonicAtomDecay ( )
virtual

Definition at line 86 of file G4MuonicAtomDecay.cc.

87{}

Member Function Documentation

◆ AtRestDoIt()

virtual G4VParticleChange * G4MuonicAtomDecay::AtRestDoIt ( const G4Track theTrack,
const G4Step theStep 
)
inlinevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 86 of file G4MuonicAtomDecay.hh.

88 {return DecayIt(theTrack, theStep);}

◆ AtRestGetPhysicalInteractionLength()

G4double G4MuonicAtomDecay::AtRestGetPhysicalInteractionLength ( const G4Track aTrack,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 113 of file G4MuonicAtomDecay.cc.

115{
117 // check if this is the beginning of tracking
120 }
122}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335

◆ GetMeanFreePath()

G4double G4MuonicAtomDecay::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 585 of file G4MuonicAtomDecay.cc.

586{
587 // based on G4Decay::GetMeanFreePath check; fixme
588
589 // get particle
590 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
591 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
592 G4double aMass = aParticle->GetMass();
593 G4double aLife = aParticleDef->GetPDGLifeTime();
594
595 // returns the mean free path in GEANT4 internal units
596 G4double pathlength;
597 G4double aCtau = c_light * aLife;
598
599 // check if the particle is stable?
600 if (aParticleDef->GetPDGStable()) {
601 pathlength = DBL_MAX;
602
603 //check if the particle has very short life time ?
604 } else if (aCtau < DBL_MIN) {
605 pathlength = DBL_MIN;
606
607 } else {
608 //calculate the mean free path
609 // by using normalized kinetic energy (= Ekin/mass)
610 G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
611 const G4double HighestValue = 20.0; //
612 if ( rKineticEnergy > HighestValue) {
613 // gamma >> 1
614 pathlength = ( rKineticEnergy + 1.0)* aCtau;
615 } else if ( rKineticEnergy < DBL_MIN ) {
616 // too slow particle
617#ifdef G4VERBOSE
618 if (GetVerboseLevel()>1) {
619 G4cout << "G4MuonicAtomDecay::GetMeanFreePath() !!particle stops!!";
620 G4cout << aParticleDef->GetParticleName() << G4endl;
621 G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
622 }
623#endif
624 pathlength = DBL_MIN;
625 } else {
626 // beta <1
627 pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
628 }
629 }
630 return pathlength;
631}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalMomentum() const
G4int GetVerboseLevel() const
G4bool GetPDGStable() const
G4double GetPDGLifeTime() const
const G4String & GetParticleName() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62

◆ GetMeanLifeTime()

G4double G4MuonicAtomDecay::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 135 of file G4MuonicAtomDecay.cc.

137{
138 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
139 G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
140 G4MuonicAtom* muatom = static_cast<G4MuonicAtom*>(aParticleDef);
141 G4double meanlife = muatom->GetPDGLifeTime();
142#ifdef G4VERBOSE
143 if (GetVerboseLevel()>1) {
144 G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
145 }
146#endif
147 return meanlife;
148}
#define ns(x)
Definition: xmltok.c:1649

Referenced by AtRestGetPhysicalInteractionLength().

◆ GetVerboseLevel()

G4int G4MuonicAtomDecay::GetVerboseLevel ( ) const
inline

Definition at line 74 of file G4MuonicAtomDecay.hh.

74{return verboseLevel;}

Referenced by GetMeanFreePath(), and GetMeanLifeTime().

◆ IsApplicable()

G4bool G4MuonicAtomDecay::IsApplicable ( const G4ParticleDefinition a)
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4MuonicAtomDecay.cc.

92{
93 return ( a.GetParticleType() == "MuonicAtom" );
94}
const G4String & GetParticleType() const

◆ PostStepDoIt()

virtual G4VParticleChange * G4MuonicAtomDecay::PostStepDoIt ( const G4Track theTrack,
const G4Step theStep 
)
inlinevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 90 of file G4MuonicAtomDecay.hh.

92 {return DecayIt(theTrack, theStep);}

◆ PostStepGetPhysicalInteractionLength()

G4double G4MuonicAtomDecay::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 126 of file G4MuonicAtomDecay.cc.

128{
130 return DBL_MAX; // this will need to be changed in future; fixme
131}

◆ ProcessDescription()

void G4MuonicAtomDecay::ProcessDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VProcess.

Definition at line 443 of file G4MuonicAtomDecay.cc.

444{
445 outFile << "MuonicAtom process where Muon decays in orbit or is captured by the nucleus." <<G4endl;
446}

◆ SetVerboseLevel()

void G4MuonicAtomDecay::SetVerboseLevel ( G4int  value)
inline

Definition at line 71 of file G4MuonicAtomDecay.hh.

71{verboseLevel = value;}

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