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

#include <G4MuonMinusCaptureAtRest.hh>

+ Inheritance diagram for G4MuonMinusCaptureAtRest:

Public Member Functions

 G4MuonMinusCaptureAtRest (const G4String &processName="muMinusCaptureAtRest", G4ProcessType aType=fHadronic)
 
virtual ~G4MuonMinusCaptureAtRest ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)
 
- Public Member Functions inherited from G4VRestProcess
 G4VRestProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestProcess (G4VRestProcess &)
 
virtual ~G4VRestProcess ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Additional Inherited Members

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

Detailed Description

Definition at line 64 of file G4MuonMinusCaptureAtRest.hh.

Constructor & Destructor Documentation

◆ G4MuonMinusCaptureAtRest()

G4MuonMinusCaptureAtRest::G4MuonMinusCaptureAtRest ( const G4String processName = "muMinusCaptureAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 61 of file G4MuonMinusCaptureAtRest.cc.

62 :
63 G4VRestProcess (processName, aType), nCascade(0), targetZ(0), targetA(0),
64 isInitialised(false)
65{
67 Cascade = new G4GHEKinematicsVector [17];
68 pSelector = new G4StopElementSelector();
69 pEMCascade = new G4MuMinusCaptureCascade();
70 theN = new G4Fancy3DNucleus();
71 theHandler = new G4ExcitationHandler();
73}
@ fHadronAtRest
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4MuonMinusCaptureAtRest()

G4MuonMinusCaptureAtRest::~G4MuonMinusCaptureAtRest ( )
virtual

Definition at line 77 of file G4MuonMinusCaptureAtRest.cc.

78{
80 delete [] Cascade;
81 delete pSelector;
82 delete pEMCascade;
83 delete theN;
84 delete theHandler;
85}
void DeRegisterExtraProcess(G4VProcess *)

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VRestProcess.

Definition at line 110 of file G4MuonMinusCaptureAtRest.cc.

112{
113 //
114 // Handles MuonMinuss at rest; a MuonMinus can either create secondaries or
115 // do nothing (in which case it should be sent back to decay-handling
116 // section
117 //
119
120 // select element and get Z,A.
121 G4Element* aEle = pSelector->GetElement(track.GetMaterial());
122 targetZ = aEle->GetZ();
123 targetA = G4double(G4int(aEle->GetN()+0.5));
124 G4int ni = 0;
125
126 G4IsotopeVector* isv = aEle->GetIsotopeVector();
127 if(isv) ni = isv->size();
128
129 if(ni == 1) {
130 targetA = G4double(aEle->GetIsotope(0)->GetN());
131 } else if(ni > 1) {
134 G4int j = -1;
135 ni--;
136 do {
137 j++;
138 y -= ab[j];
139 } while (y > 0.0 && j < ni);
140 targetA = G4double(aEle->GetIsotope(j)->GetN());
141 }
142
143 // Do the electromagnetic cascade of the muon in the nuclear field.
144 nCascade = 0;
145 targetMass = G4NucleiProperties::GetNuclearMass(targetA, targetZ);
146 nCascade = pEMCascade->DoCascade(targetZ, targetMass, Cascade);
147
148 // Decide on Decay or Capture, and doit.
149 G4double lambdac = pSelector->GetMuonCaptureRate(targetZ, targetA);
150 G4double lambdad = pSelector->GetMuonDecayRate(targetZ, targetA);
151 G4double lambda = lambdac + lambdad;
152
153 // === Throw for capture time.
154
155 G4double tDelay = -std::log(G4UniformRand()) / lambda;
156
157 G4ReactionProductVector * captureResult = 0;
158 G4int nEmSecondaries = nCascade;
159 G4int nSecondaries = nCascade;
160 /*
161 G4cout << "lambda= " << lambda << " lambdac= " << lambdac
162 << " nem= " << nEmSecondaries << G4endl;
163 */
164 if( G4UniformRand()*lambda > lambdac)
165 pEMCascade->DoBoundMuonMinusDecay(targetZ, &nEmSecondaries, Cascade);
166 else
167 captureResult = DoMuCapture();
168
169 // fill the final state
170 if(captureResult) nSecondaries += captureResult->size();
171 else nSecondaries = nEmSecondaries;
172 //G4cout << " nsec= " << nSecondaries << " nem= " << nEmSecondaries << G4endl;
173
175
176 G4double globalTime = track.GetGlobalTime();
178 // Store nuclear cascade
179 if(captureResult) {
180 G4int n = captureResult->size();
181 for ( G4int isec = 0; isec < n; isec++ ) {
182 G4ReactionProduct* aParticle = captureResult->operator[](isec);
183 G4DynamicParticle * aNewParticle = new G4DynamicParticle();
184 aNewParticle->SetDefinition( aParticle->GetDefinition() );
185 G4LorentzVector itV(aParticle->GetTotalEnergy(), aParticle->GetMomentum());
186 aNewParticle->SetMomentum(itV.vect());
187 G4double localtime = globalTime + tDelay + aParticle->GetTOF();
188 G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position);
189 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
190 aParticleChange.AddSecondary( aNewTrack );
191 delete aParticle;
192 }
193 delete captureResult;
194 }
195
196 // Store electromagnetic cascade
197
198 if(nEmSecondaries > 0) {
199
200 for ( G4int isec = 0; isec < nEmSecondaries; isec++ ) {
201 G4ParticleDefinition* pd = Cascade[isec].GetParticleDef();
202 G4double localtime = globalTime;
203 if(isec >= nCascade) localtime += tDelay;
204 if(pd) {
205 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
206 aNewParticle->SetDefinition( pd );
207 aNewParticle->SetMomentum( Cascade[isec].GetMomentum() );
208
209 G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position );
210 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
211 aParticleChange.AddSecondary( aNewTrack );
212 }
213 }
214 }
215
218
219 return &aParticleChange;
220}
std::vector< G4Isotope * > G4IsotopeVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:131
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4double GetN() const
Definition: G4Element.hh:134
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
G4ParticleDefinition * GetParticleDef()
G4int GetN() const
Definition: G4Isotope.hh:94
void DoBoundMuonMinusDecay(G4double Z, G4int *nCascade, G4GHEKinematicsVector *Cascade)
G4int DoCascade(const G4double Z, const G4double A, G4GHEKinematicsVector *Cascade)
static G4double GetNuclearMass(const G4double A, const G4double Z)
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTOF() const
G4ParticleDefinition * GetDefinition() const
G4Element * GetElement(const G4Material *aMaterial)
G4double GetMuonDecayRate(G4double Z, G4double A)
G4double GetMuonCaptureRate(G4double Z, G4double A)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289

◆ BuildPhysicsTable()

void G4MuonMinusCaptureAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 103 of file G4MuonMinusCaptureAtRest.cc.

104{
106}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

G4double G4MuonMinusCaptureAtRest::GetMeanLifeTime ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VRestProcess.

Definition at line 82 of file G4MuonMinusCaptureAtRest.hh.

83 {return 0;};

◆ IsApplicable()

G4bool G4MuonMinusCaptureAtRest::IsApplicable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 89 of file G4MuonMinusCaptureAtRest.cc.

90{
91 return ( &p == G4MuonMinus::MuonMinus() );
92}
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100

◆ PreparePhysicsTable()

void G4MuonMinusCaptureAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 96 of file G4MuonMinusCaptureAtRest.cc.

97{
99}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)

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