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

#include <G4Decay.hh>

+ Inheritance diagram for G4Decay:

Public Member Functions

 G4Decay (const G4String &processName="Decay")
 
virtual ~G4Decay ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual void StartTracking (G4Track *) override
 
virtual void EndTracking () override
 
void SetExtDecayer (G4VExtDecayer *)
 
const G4VExtDecayerGetExtDecayer () const
 
G4double GetRemainderLifeTime () const
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
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
 
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 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 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
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual G4VParticleChangeDecayIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void DaughterPolarization (const G4Track &aTrack, G4DecayProducts *products)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition) override
 
- Protected Member Functions inherited from G4VRestDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4int verboseLevel
 
const G4double HighestValue
 
G4double fRemainderLifeTime
 
G4ParticleChangeForDecay fParticleChangeForDecay
 
G4VExtDecayerpExtDecayer
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 62 of file G4Decay.hh.

Constructor & Destructor Documentation

◆ G4Decay()

G4Decay::G4Decay ( const G4String & processName = "Decay")

Definition at line 62 of file G4Decay.cc.

63 :G4VRestDiscreteProcess(processName, fDecay),
64 verboseLevel(1),
65 HighestValue(20.0),
67 pExtDecayer(nullptr)
68{
69 // set Process Sub Type
70 SetProcessSubType(static_cast<int>(DECAY));
71
72#ifdef G4VERBOSE
73 if (GetVerboseLevel()>1) {
74 G4cout << "G4Decay constructor " << " Name:" << processName << G4endl;
75 }
76#endif
77
79}
@ fDecay
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VExtDecayer * pExtDecayer
Definition G4Decay.hh:179
G4double fRemainderLifeTime
Definition G4Decay.hh:173
G4ParticleChangeForDecay fParticleChangeForDecay
Definition G4Decay.hh:176
const G4double HighestValue
Definition G4Decay.hh:170
G4int verboseLevel
Definition G4Decay.hh:162
G4int GetVerboseLevel() const
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange

◆ ~G4Decay()

G4Decay::~G4Decay ( )
virtual

Definition at line 81 of file G4Decay.cc.

82{
83 if (pExtDecayer != nullptr) {
84 delete pExtDecayer;
85 }
86}

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4Decay::AtRestDoIt ( const G4Track & aTrack,
const G4Step & aStep )
inlineoverridevirtual

Reimplemented from G4VRestDiscreteProcess.

Reimplemented in G4DecayWithSpin.

Definition at line 183 of file G4Decay.hh.

187{
188 return DecayIt(aTrack, aStep);
189}
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition G4Decay.cc:180

◆ AtRestGetPhysicalInteractionLength()

G4double G4Decay::AtRestGetPhysicalInteractionLength ( const G4Track & track,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 477 of file G4Decay.cc.

481{
482 // condition is set to "Not Forced"
484
486 if (pTime >= 0.) {
487 fRemainderLifeTime = pTime - track.GetProperTime();
489 } else {
492 }
493 return fRemainderLifeTime;
494}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition G4Types.hh:83
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition) override
Definition G4Decay.cc:100
G4double GetPreAssignedDecayProperTime() const
G4double GetProperTime() const
const G4DynamicParticle * GetDynamicParticle() const
G4double theNumberOfInteractionLengthLeft
#define DBL_MIN
Definition templates.hh:54

◆ BuildPhysicsTable()

void G4Decay::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 175 of file G4Decay.cc.

176{
177 return;
178}

◆ DaughterPolarization()

void G4Decay::DaughterPolarization ( const G4Track & aTrack,
G4DecayProducts * products )
protectedvirtual

Reimplemented in G4PionDecayMakeSpin.

Definition at line 386 of file G4Decay.cc.

387{
388 // empty implementation
389}

Referenced by DecayIt().

◆ DecayIt()

G4VParticleChange * G4Decay::DecayIt ( const G4Track & aTrack,
const G4Step & aStep )
protectedvirtual

Definition at line 180 of file G4Decay.cc.

181{
182 // The DecayIt() method returns by pointer a particle-change object.
183 // Units are expressed in GEANT4 internal units.
184
185 // Initialize ParticleChange
186 // all members of G4VParticleChange are set to equal to
187 // corresponding member in G4Track
189
190 // get particle
191 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
192 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
193
194 // check if the particle is stable
195 if (aParticleDef->GetPDGStable()) return &fParticleChangeForDecay ;
196
197
198 //check if thePreAssignedDecayProducts exists
199 const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
200 G4bool isPreAssigned = (o_products != nullptr);
201 G4DecayProducts* products = nullptr;
202
203 // decay table
204 G4DecayTable *decaytable = aParticleDef->GetDecayTable();
205
206 // check if external decayer exists
207 G4bool isExtDecayer = (decaytable == nullptr) && (pExtDecayer != nullptr);
208
209 // Error due to NO Decay Table
210 if ( (decaytable == nullptr) && !isExtDecayer && !isPreAssigned ){
211 if (GetVerboseLevel()>0) {
212 G4cout << "G4Decay::DoIt : decay table not defined for ";
213 G4cout << aParticle->GetDefinition()->GetParticleName()<< G4endl;
214 }
216 ed << "For " << aParticle->GetDefinition()->GetParticleName()
217 << " decay probability exist but decay table is not defined "
218 << "- the particle will be killed;\n"
219 << " isExtDecayer: " << isExtDecayer
220 << "; isPreAssigned: " << isPreAssigned;
221 G4Exception( "G4Decay::DecayIt ",
222 "DECAY101",JustWarning, ed);
223
225 // Kill the parent particle
228
231 }
232
233 if (isPreAssigned) {
234 // copy decay products
235 products = new G4DecayProducts(*o_products);
236 } else if ( isExtDecayer ) {
237 // decay according to external decayer
238 products = pExtDecayer->ImportDecayProducts(aTrack);
239 } else {
240 // Decay according to decay table.
241 // Keep trying to choose a candidate decay channel if the dynamic mass
242 // of the decaying particle is below the sum of the PDG masses of the
243 // candidate daughter particles.
244 // This is needed because the decay table used in Geant4 is based on
245 // the assumption of nominal PDG masses, but a wide resonance can have
246 // a dynamic masses well below its nominal PDG masses, and therefore
247 // some of its decay channels can be below the kinematical threshold.
248 // Note that, for simplicity, we ignore here the possibility that
249 // one or more of the candidate daughter particles can be, in turn,
250 // wide resonance. However, if this is the case, and the channel is
251 // accepted, then the masses of the resonance daughter particles will
252 // be sampled by taking into account their widths.
253 G4VDecayChannel* decaychannel = nullptr;
254 G4double massParent = aParticle->GetMass();
255 decaychannel = decaytable->SelectADecayChannel(massParent);
256 if ( decaychannel == nullptr) {
257 // decay channel not found
259 ed << "Can not determine decay channel for "
260 << aParticleDef->GetParticleName() << G4endl
261 << " mass of dynamic particle: "
262 << massParent/GeV << " (GEV)" << G4endl
263 << " dacay table has " << decaytable->entries()
264 << " entries" << G4endl;
265 G4double checkedmass=massParent;
266 if (massParent < 0.) {
267 checkedmass=aParticleDef->GetPDGMass();
268 ed << "Using PDG mass ("<<checkedmass/GeV
269 << "(GeV)) in IsOKWithParentMass" << G4endl;
270 }
271 for (G4int ic =0;ic <decaytable->entries();++ic) {
272 G4VDecayChannel * dc= decaytable->GetDecayChannel(ic);
273 ed << ic << ": BR " << dc->GetBR() << ", IsOK? "
274 << dc->IsOKWithParentMass(checkedmass)
275 << ", --> ";
276 G4int ndaughters=dc->GetNumberOfDaughters();
277 for (G4int id=0;id<ndaughters;++id) {
278 if (id>0) ed << " + "; // seperator, except for first
279 ed << dc->GetDaughterName(id);
280 }
281 ed << G4endl;
282 }
283 G4Exception("G4Decay::DoIt", "DECAY003", FatalException,ed);
284 } else {
285 // execute DecayIt()
286#ifdef G4VERBOSE
287 G4int temp = decaychannel->GetVerboseLevel();
288 if (GetVerboseLevel()>1) {
289 G4cout << "G4Decay::DoIt : selected decay channel addr:"
290 << decaychannel <<G4endl;
291 decaychannel->SetVerboseLevel(GetVerboseLevel());
292 }
293#endif
294 products = decaychannel->DecayIt(aParticle->GetMass());
295#ifdef G4VERBOSE
296 if (GetVerboseLevel()>1) {
297 decaychannel->SetVerboseLevel(temp);
298 }
299#endif
300#ifdef G4VERBOSE
301 if (GetVerboseLevel()>2) {
302 if (! products->IsChecked() ) products->DumpInfo();
303 }
304#endif
305 }
306 }
307
308 // get parent particle information ...................................
309 G4double ParentEnergy = aParticle->GetTotalEnergy();
310 G4double ParentMass = aParticle->GetMass();
311 if (ParentEnergy < ParentMass) {
313 ed << "Total Energy is less than its mass - increased the energy"
314 << "\n Particle: " << aParticle->GetDefinition()->GetParticleName()
315 << "\n Energy:" << ParentEnergy/MeV << "[MeV]"
316 << "\n Mass:" << ParentMass/MeV << "[MeV]";
317 G4Exception( "G4Decay::DecayIt ",
318 "DECAY102",JustWarning, ed);
319 ParentEnergy = ParentMass;
320 }
321
322 G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
323
324 //boost all decay products to laboratory frame
325 G4double energyDeposit = 0.0;
326 G4double finalGlobalTime = aTrack.GetGlobalTime();
327 G4double finalLocalTime = aTrack.GetLocalTime();
328 if (aTrack.GetTrackStatus() == fStopButAlive ){
329 // AtRest case
330 finalGlobalTime += fRemainderLifeTime;
331 finalLocalTime += fRemainderLifeTime;
332 energyDeposit += aParticle->GetKineticEnergy();
333 if (isPreAssigned) products->Boost( ParentEnergy, ParentDirection);
334 } else {
335 // PostStep case
336 if (!isExtDecayer) products->Boost( ParentEnergy, ParentDirection);
337 }
338 // set polarization for daughter particles
339 DaughterPolarization(aTrack, products);
340
341
342 //add products in fParticleChangeForDecay
343 G4int numberOfSecondaries = products->entries();
345#ifdef G4VERBOSE
346 if (GetVerboseLevel()>1) {
347 G4cout << "G4Decay::DoIt : Decay vertex :";
348 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
349 G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
350 G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
351 G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
352 G4cout << G4endl;
353 G4cout << "G4Decay::DoIt : decay products in Lab. Frame" << G4endl;
354 products->DumpInfo();
355 }
356#endif
357 G4int index;
358 G4ThreeVector currentPosition;
359 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
360 for (index=0; index < numberOfSecondaries; index++){
361 // get current position of the track
362 currentPosition = aTrack.GetPosition();
363 // create a new track object
364 G4Track* secondary = new G4Track( products->PopProducts(),
365 finalGlobalTime ,
366 currentPosition );
367 // switch on good for tracking flag
368 secondary->SetGoodForTrackingFlag();
369 secondary->SetTouchableHandle(thand);
370 // add the secondary track in the List
372 }
373 delete products;
374
375 // Kill the parent particle
379
380 // Clear NumberOfInteractionLengthLeft
382
384}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ fStopAndKill
@ fStopButAlive
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
G4bool IsChecked() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * GetDecayChannel(G4int index) const
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
G4int entries() const
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition G4Decay.cc:386
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void Initialize(const G4Track &) final
G4bool GetPDGStable() const
G4DecayTable * GetDecayTable() const
const G4String & GetParticleName() const
G4TrackStatus GetTrackStatus() const
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4double GetLocalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void SetGoodForTrackingFlag(G4bool value=true)
G4double GetBR() const
void SetVerboseLevel(G4int value)
G4int GetVerboseLevel() const
G4int GetNumberOfDaughters() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
virtual G4bool IsOKWithParentMass(G4double parentMass)
const G4String & GetDaughterName(G4int anIndex) const
virtual G4DecayProducts * ImportDecayProducts(const G4Track &aTrack)=0
void ProposeTrackStatus(G4TrackStatus status)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
void ClearNumberOfInteractionLengthLeft()

Referenced by AtRestDoIt(), G4DecayWithSpin::AtRestDoIt(), PostStepDoIt(), and G4DecayWithSpin::PostStepDoIt().

◆ EndTracking()

void G4Decay::EndTracking ( )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 401 of file G4Decay.cc.

402{
403 // Clear NumberOfInteractionLengthLeft
405
407}
G4double currentInteractionLength

◆ GetExtDecayer()

const G4VExtDecayer * G4Decay::GetExtDecayer ( ) const
inline

Definition at line 193 of file G4Decay.hh.

194{
195 return pExtDecayer;
196}

◆ GetMeanFreePath()

G4double G4Decay::GetMeanFreePath ( const G4Track & aTrack,
G4double previousStepSize,
G4ForceCondition * condition )
overrideprotectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 129 of file G4Decay.cc.

130{
131 // get particle
132 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
133 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
134 G4double aMass = aParticle->GetMass();
135 G4double aLife = aParticleDef->GetPDGLifeTime();
136
137
138 // returns the mean free path in GEANT4 internal units
139 G4double pathlength;
140 G4double aCtau = c_light * aLife;
141
142 // check if the particle is stable?
143 if (aParticleDef->GetPDGStable()) {
144 pathlength = DBL_MAX;
145
146 //check if the particle has very short life time ?
147 } else if (aCtau < DBL_MIN) {
148 pathlength = DBL_MIN;
149
150 } else {
151 //calculate the mean free path
152 // by using normalized kinetic energy (= Ekin/mass)
153 G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
154 if ( rKineticEnergy > HighestValue) {
155 // gamma >> 1
156 pathlength = ( rKineticEnergy + 1.0)* aCtau;
157 } else if ( rKineticEnergy < DBL_MIN ) {
158 // too slow particle
159#ifdef G4VERBOSE
160 if (GetVerboseLevel()>1) {
161 G4cout << "G4Decay::GetMeanFreePath() !!particle stops!!";
162 G4cout << aParticleDef->GetParticleName() << G4endl;
163 G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
164 }
165#endif
166 pathlength = DBL_MIN;
167 } else {
168 // beta <1
169 pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
170 }
171 }
172 return pathlength;
173}
G4double GetTotalMomentum() const
G4double GetPDGLifeTime() const
#define DBL_MAX
Definition templates.hh:62

Referenced by PostStepGetPhysicalInteractionLength().

◆ GetMeanLifeTime()

G4double G4Decay::GetMeanLifeTime ( const G4Track & aTrack,
G4ForceCondition * condition )
overrideprotectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 100 of file G4Decay.cc.

102{
103 // returns the mean free path in GEANT4 internal units
104 G4double meanlife;
105
106 // get particle
107 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
108 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
109 G4double aLife = aParticleDef->GetPDGLifeTime();
110
111 // check if the particle is stable?
112 if (aParticleDef->GetPDGStable()) {
113 //1000000 times the life time of the universe
114 meanlife = 1e24 * s;
115
116 } else {
117 meanlife = aLife;
118 }
119
120#ifdef G4VERBOSE
121 if (GetVerboseLevel()>1) {
122 G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
123 }
124#endif
125
126 return meanlife;
127}

Referenced by AtRestGetPhysicalInteractionLength().

◆ GetRemainderLifeTime()

G4double G4Decay::GetRemainderLifeTime ( ) const
inline

Definition at line 199 of file G4Decay.hh.

200{
201 return fRemainderLifeTime;
202}

◆ IsApplicable()

G4bool G4Decay::IsApplicable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 88 of file G4Decay.cc.

89{
90 // check if the particle is stable?
91 if (aParticleType.GetPDGLifeTime() <0.0) {
92 return false;
93 } else if (aParticleType.GetPDGMass() <= 0.0*MeV) {
94 return false;
95 } else {
96 return true;
97 }
98}

Referenced by LBE::ConstructGeneral(), and G4DecayPhysics::ConstructProcess().

◆ PostStepDoIt()

G4VParticleChange * G4Decay::PostStepDoIt ( const G4Track & aTrack,
const G4Step & aStep )
overridevirtual

Reimplemented from G4VRestDiscreteProcess.

Reimplemented in G4DecayWithSpin.

Definition at line 507 of file G4Decay.cc.

511{
512 if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
513 (aTrack.GetTrackStatus() == fStopAndKill ) ){
516 } else {
517 return DecayIt(aTrack, aStep);
518 }
519}

◆ PostStepGetPhysicalInteractionLength()

G4double G4Decay::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 410 of file G4Decay.cc.

415{
416 // condition is set to "Not Forced"
418
419 // pre-assigned Decay time
422
423 if (pTime < 0.) {
424 // normal case
425 if ( previousStepSize > 0.0){
426 // subtract NumberOfInteractionLengthLeft
430 }
432 }
433 // get mean free path
434 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
435
436#ifdef G4VERBOSE
437 if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
438 G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
439 track.GetDynamicParticle()->DumpInfo();
440 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
441 G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
442 }
443#endif
444
445 G4double value;
448 //fRemainderLifeTime = theNumberOfInteractionLengthLeft*aLife;
449 } else {
450 value = DBL_MAX;
451 }
452
453 return value;
454
455 } else {
456 //pre-assigned Decay time case
457 // reminder proper time
458 fRemainderLifeTime = pTime - track.GetProperTime();
459 if (fRemainderLifeTime <= 0.0) fRemainderLifeTime = 0.0;
460
461 G4double rvalue=0.0;
462 // use pre-assigned Decay time to determine PIL
463 if (aLife>0.0) {
464 // ordinary particle
465 rvalue = (fRemainderLifeTime/aLife)*GetMeanFreePath(track, previousStepSize, condition);
466 } else {
467 // shortlived particle
468 rvalue = c_light * fRemainderLifeTime;
469 // by using normalized kinetic energy (= Ekin/mass)
470 G4double aMass = track.GetDynamicParticle()->GetMass();
471 rvalue *= track.GetDynamicParticle()->GetTotalMomentum()/aMass;
472 }
473 return rvalue;
474 }
475}
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
Definition G4Decay.cc:129
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
G4Material * GetMaterial() const
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)

◆ ProcessDescription()

void G4Decay::ProcessDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4DecayWithSpin, and G4PionDecayMakeSpin.

Definition at line 521 of file G4Decay.cc.

522{
523 outFile << GetProcessName() << ": Decay of particles. \n"
524 << "kinematics of daughters are dertermined by DecayChannels "
525 << " or by PreAssignedDecayProducts\n";
526}
const G4String & GetProcessName() const

◆ SetExtDecayer()

void G4Decay::SetExtDecayer ( G4VExtDecayer * val)

Definition at line 497 of file G4Decay.cc.

498{
499 pExtDecayer = val;
500
501 // set Process Sub Type
502 if ( pExtDecayer !=0 ) {
503 SetProcessSubType(static_cast<int>(DECAY_External));
504 }
505}
@ DECAY_External

◆ StartTracking()

void G4Decay::StartTracking ( G4Track * )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 393 of file G4Decay.cc.

394{
397
398 fRemainderLifeTime = -1.0;
399}
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80

Member Data Documentation

◆ fParticleChangeForDecay

G4ParticleChangeForDecay G4Decay::fParticleChangeForDecay
protected

Definition at line 176 of file G4Decay.hh.

Referenced by DecayIt(), G4Decay(), PostStepDoIt(), and G4DecayWithSpin::PostStepDoIt().

◆ fRemainderLifeTime

◆ HighestValue

const G4double G4Decay::HighestValue
protected

Definition at line 170 of file G4Decay.hh.

Referenced by GetMeanFreePath().

◆ pExtDecayer

G4VExtDecayer* G4Decay::pExtDecayer
protected

Definition at line 179 of file G4Decay.hh.

Referenced by DecayIt(), GetExtDecayer(), SetExtDecayer(), and ~G4Decay().

◆ verboseLevel

G4int G4Decay::verboseLevel
protected

Definition at line 162 of file G4Decay.hh.

Referenced by PostStepGetPhysicalInteractionLength().


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