Geant4 11.4.0
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 G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Protected Attributes

G4int verboseLevel
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),
66 pExtDecayer(nullptr)
67{
68 // set Process Sub Type
69 SetProcessSubType(static_cast<int>(DECAY));
70
71#ifdef G4VERBOSE
72 if (GetVerboseLevel()>1) {
73 G4cout << "G4Decay constructor " << " Name:" << processName << G4endl;
74 }
75#endif
76
78}
@ fDecay
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VExtDecayer * pExtDecayer
Definition G4Decay.hh:177
G4double fRemainderLifeTime
Definition G4Decay.hh:171
G4ParticleChangeForDecay fParticleChangeForDecay
Definition G4Decay.hh:174
G4int verboseLevel
Definition G4Decay.hh:162
G4int GetVerboseLevel() const
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
G4VRestDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)

Referenced by G4DecayWithSpin::G4DecayWithSpin(), and G4PionDecayMakeSpin::G4PionDecayMakeSpin().

◆ ~G4Decay()

G4Decay::~G4Decay ( )
virtual

Definition at line 80 of file G4Decay.cc.

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

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Reimplemented in G4DecayWithSpin.

Definition at line 181 of file G4Decay.hh.

185{
186 return DecayIt(aTrack, aStep);
187}
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition G4Decay.cc:176

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 473 of file G4Decay.cc.

477{
478 // condition is set to "Not Forced"
480
482 if (pTime >= 0.) {
483 fRemainderLifeTime = pTime - track.GetProperTime();
485 } else {
488 }
489 return fRemainderLifeTime;
490}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition G4Types.hh:83
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition) override
Definition G4Decay.cc:99
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 171 of file G4Decay.cc.

172{
173 return;
174}

◆ DaughterPolarization()

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

Reimplemented in G4PionDecayMakeSpin.

Definition at line 382 of file G4Decay.cc.

383{
384 // empty implementation
385}

Referenced by DecayIt().

◆ DecayIt()

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

Definition at line 176 of file G4Decay.cc.

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

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

◆ EndTracking()

void G4Decay::EndTracking ( )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 397 of file G4Decay.cc.

398{
399 // Clear NumberOfInteractionLengthLeft
401
403}
G4double currentInteractionLength

◆ GetExtDecayer()

const G4VExtDecayer * G4Decay::GetExtDecayer ( ) const
inline

Definition at line 191 of file G4Decay.hh.

192{
193 return pExtDecayer;
194}

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 128 of file G4Decay.cc.

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

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

Referenced by AtRestGetPhysicalInteractionLength().

◆ GetRemainderLifeTime()

G4double G4Decay::GetRemainderLifeTime ( ) const
inline

Definition at line 197 of file G4Decay.hh.

198{
199 return fRemainderLifeTime;
200}

◆ IsApplicable()

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

Reimplemented from G4VProcess.

Definition at line 87 of file G4Decay.cc.

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

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 503 of file G4Decay.cc.

507{
508 if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
509 (aTrack.GetTrackStatus() == fStopAndKill ) ){
510 fParticleChangeForDecay.Initialize(aTrack);
512 } else {
513 return DecayIt(aTrack, aStep);
514 }
515}

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 406 of file G4Decay.cc.

411{
412 // condition is set to "Not Forced"
414
415 // pre-assigned Decay time
418
419 if (pTime < 0.) {
420 // normal case
421 if ( previousStepSize > 0.0){
422 // subtract NumberOfInteractionLengthLeft
426 }
428 }
429 // get mean free path
430 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
431
432#ifdef G4VERBOSE
433 if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
434 G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
435 track.GetDynamicParticle()->DumpInfo();
436 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
437 G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
438 }
439#endif
440
441 G4double value;
444 //fRemainderLifeTime = theNumberOfInteractionLengthLeft*aLife;
445 } else {
446 value = DBL_MAX;
447 }
448
449 return value;
450
451 } else {
452 //pre-assigned Decay time case
453 // reminder proper time
454 fRemainderLifeTime = pTime - track.GetProperTime();
455 if (fRemainderLifeTime <= 0.0) fRemainderLifeTime = 0.0;
456
457 G4double rvalue=0.0;
458 // use pre-assigned Decay time to determine PIL
459 if (aLife>0.0) {
460 // ordinary particle
461 rvalue = (fRemainderLifeTime/aLife)*GetMeanFreePath(track, previousStepSize, condition);
462 } else {
463 // shortlived particle
464 rvalue = c_light * fRemainderLifeTime;
465 // by using normalized kinetic energy (= Ekin/mass)
466 G4double aMass = track.GetDynamicParticle()->GetMass();
467 rvalue *= track.GetDynamicParticle()->GetTotalMomentum()/aMass;
468 }
469 return rvalue;
470 }
471}
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
Definition G4Decay.cc:128
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 517 of file G4Decay.cc.

518{
519 outFile << GetProcessName() << ": Decay of particles. \n"
520 << "kinematics of daughters are dertermined by DecayChannels "
521 << " or by PreAssignedDecayProducts\n";
522}
const G4String & GetProcessName() const

◆ SetExtDecayer()

void G4Decay::SetExtDecayer ( G4VExtDecayer * val)

Definition at line 493 of file G4Decay.cc.

494{
495 pExtDecayer = val;
496
497 // set Process Sub Type
498 if ( pExtDecayer !=0 ) {
499 SetProcessSubType(static_cast<int>(DECAY_External));
500 }
501}
@ DECAY_External

◆ StartTracking()

void G4Decay::StartTracking ( G4Track * )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 389 of file G4Decay.cc.

390{
393
394 fRemainderLifeTime = -1.0;
395}
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80

Member Data Documentation

◆ fParticleChangeForDecay

G4ParticleChangeForDecay G4Decay::fParticleChangeForDecay
protected

Definition at line 174 of file G4Decay.hh.

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

◆ fRemainderLifeTime

◆ pExtDecayer

G4VExtDecayer* G4Decay::pExtDecayer
protected

Definition at line 177 of file G4Decay.hh.

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

◆ verboseLevel

G4int G4Decay::verboseLevel
protected

Definition at line 162 of file G4Decay.hh.

Referenced by G4Decay(), and PostStepGetPhysicalInteractionLength().


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