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

#include <G4MuonMinusAtomicCapture.hh>

+ Inheritance diagram for G4MuonMinusAtomicCapture:

Public Member Functions

 G4MuonMinusAtomicCapture (const G4String &name="muMinusAtomicCaptureAtRest")
 
 ~G4MuonMinusAtomicCapture ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
void ProcessDescription (std::ostream &outFile) const
 
void SetElementSelector (G4ElementSelector *ptr)
 
void SetEmCascade (G4HadronicInteraction *ptr)
 
- Public Member Functions inherited from G4VRestProcess
 G4VRestProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VRestProcess (G4VRestProcess &)
 
virtual ~G4VRestProcess ()
 
G4VRestProcessoperator= (const G4VRestProcess &)=delete
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 

Protected Member Functions

G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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 60 of file G4MuonMinusAtomicCapture.hh.

Constructor & Destructor Documentation

◆ G4MuonMinusAtomicCapture()

G4MuonMinusAtomicCapture::G4MuonMinusAtomicCapture ( const G4String name = "muMinusAtomicCaptureAtRest")
explicit

Definition at line 63 of file G4MuonMinusAtomicCapture.cc.

65 fElementSelector(new G4ElementSelector()),
66 fEmCascade(new G4EmCaptureCascade()), // Owned by InteractionRegistry
67 theTotalResult(new G4ParticleChange()),
68 result(nullptr)
69{
72}
@ fMuAtomicCapture
@ fHadronic
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4MuonMinusAtomicCapture()

G4MuonMinusAtomicCapture::~G4MuonMinusAtomicCapture ( )

Definition at line 76 of file G4MuonMinusAtomicCapture.cc.

77{
79 delete theTotalResult;
80}
void DeRegisterExtraProcess(G4VProcess *)

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4MuonMinusAtomicCapture::AtRestDoIt ( const G4Track track,
const G4Step  
)
virtual

Reimplemented from G4VRestProcess.

Definition at line 114 of file G4MuonMinusAtomicCapture.cc.

116{
117 // if primary is not Alive then do nothing (how?)
118 theTotalResult->Initialize(track);
119
120 G4Nucleus* nucleus = &targetNucleus;
121 // the call below actually sets the nucleus params;
122 // G4Nucleus targetNucleus; is a member of G4HadronicProcess
123 // G4Element* elm =
124 fElementSelector->SelectZandA(track, nucleus);
125
126 thePro.Initialise(track); // thePro was G4HadProjectile from G4HadronicProcess
127
128 // save track time an dstart capture from zero time
129 thePro.SetGlobalTime(0.0);
130 G4double time0 = track.GetGlobalTime();
131
132 // Do the electromagnetic cascade in the nuclear field.
133 // EM cascade should keep G4HadFinalState object,
134 // because it will not be deleted at the end of this method
135 //
136 result = fEmCascade->ApplyYourself(thePro, *nucleus);
137 G4double ebound = result->GetLocalEnergyDeposit(); // may need to carry this over; review
138 G4double edep = 0.0;
139 G4int nSecondaries = (G4int)result->GetNumberOfSecondaries();
140 thePro.SetBoundEnergy(ebound);
141
142 // creating the muonic atom
143 ++nSecondaries;
144
146 G4ParticleDefinition* muonicAtom = itp->GetMuonicAtom(nucleus->GetZ_asInt(),
147 nucleus->GetA_asInt());
148
149 G4DynamicParticle* dp = new G4DynamicParticle(muonicAtom,G4RandomDirection(),0.);
150 G4HadSecondary hadSec(dp);
151 hadSec.SetTime(time0);
152 result->AddSecondary(hadSec);
153
154 // Fill results
155 //
156 theTotalResult->ProposeTrackStatus(fStopAndKill);
157 theTotalResult->ProposeLocalEnergyDeposit(edep);
158 theTotalResult->SetNumberOfSecondaries(nSecondaries);
159 G4double w = track.GetWeight();
160 theTotalResult->ProposeWeight(w);
161
162#ifdef G4VERBOSE
163 if (GetVerboseLevel() > 1) {
164 G4cout << __func__
165 << " nSecondaries "
166 << nSecondaries
167 << G4endl;
168 }
169#endif
170
171 for(G4int i=0; i<nSecondaries; ++i) {
172 G4HadSecondary* sec = result->GetSecondary(i);
173
174 // add track global time to the reaction time
175 G4double time = sec->GetTime();
176 if(time < 0.0) { time = 0.0; }
177 time += time0;
178
179#ifdef G4VERBOSE
180 if (GetVerboseLevel() > 1) {
181 G4cout << __func__
182 << " "
183 << i
184 << " Resulting secondary "
185 << sec->GetParticle()->GetPDGcode()
186 << " "
188 << G4endl;
189 }
190#endif
191
192 // create secondary track
193 G4Track* t = new G4Track(sec->GetParticle(),
194 time,
195 track.GetPosition());
196 t->SetWeight(w*sec->GetWeight());
197
199 theTotalResult->AddSecondary(t);
200 }
201 result->Clear();
202
203 // fixme: needs to be done at the MuonicAtom level
204 // if (epReportLevel != 0) { // G4HadronicProcess::
205 // CheckEnergyMomentumConservation(track, *nucleus);
206 // }
207 return theTotalResult;
208}
G4ThreeVector G4RandomDirection()
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4int GetPDGcode() const
virtual const G4Element * SelectZandA(const G4Track &track, G4Nucleus *)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
G4double GetLocalEnergyDeposit() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void Initialise(const G4Track &aT)
void SetGlobalTime(G4double t)
void SetBoundEnergy(G4double e)
G4DynamicParticle * GetParticle()
G4double GetWeight() const
G4double GetTime() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ParticleDefinition * GetMuonicAtom(G4Ions const *)
Definition: G4IonTable.cc:2167
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
const G4String & GetParticleName() const
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4int GetVerboseLevel() const
Definition: G4VProcess.hh:422

◆ AtRestGetPhysicalInteractionLength()

G4double G4MuonMinusAtomicCapture::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VRestProcess.

Definition at line 105 of file G4MuonMinusAtomicCapture.cc.

107{
109 return 0.0;
110}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ BuildPhysicsTable()

void G4MuonMinusAtomicCapture::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 98 of file G4MuonMinusAtomicCapture.cc.

99{
101}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

G4double G4MuonMinusAtomicCapture::GetMeanLifeTime ( const G4Track ,
G4ForceCondition  
)
inlineprotectedvirtual

Implements G4VRestProcess.

Definition at line 91 of file G4MuonMinusAtomicCapture.hh.

92 { return -1.0; }

◆ IsApplicable()

G4bool G4MuonMinusAtomicCapture::IsApplicable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 84 of file G4MuonMinusAtomicCapture.cc.

85{
86 return (&p == G4MuonMinus::MuonMinus());
87}
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99

◆ PreparePhysicsTable()

void G4MuonMinusAtomicCapture::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4MuonMinusAtomicCapture.cc.

92{
94}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 212 of file G4MuonMinusAtomicCapture.cc.

213{
214 outFile << "Stopping of mu- using default element selector, EM cascade"
215 << "G4MuonicAtom is created\n";
216}

◆ SetElementSelector()

void G4MuonMinusAtomicCapture::SetElementSelector ( G4ElementSelector ptr)
inline

◆ SetEmCascade()

void G4MuonMinusAtomicCapture::SetEmCascade ( G4HadronicInteraction ptr)
inline

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