Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Decay Class Reference

#include <G4Decay.hh>

Inheritance diagram for G4Decay:
Collaboration diagram for G4Decay:

Public Member Functions

 G4Decay (const G4String &processName="Decay")
 
virtual ~G4Decay ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
void SetExtDecayer (G4VExtDecayer *)
 
const G4VExtDecayerGetExtDecayer () const
 
G4double GetRemainderLifeTime () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

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

Additional Inherited Members

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

Detailed Description

Definition at line 63 of file G4Decay.hh.

Constructor & Destructor Documentation

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

Definition at line 63 of file G4Decay.cc.

64  :G4VRestDiscreteProcess(processName, fDecay),
65  verboseLevel(1),
66  HighestValue(20.0),
67  fRemainderLifeTime(-1.0),
68  pExtDecayer(0)
69 {
70  // set Process Sub Type
71  SetProcessSubType(static_cast<int>(DECAY));
72 
73 #ifdef G4VERBOSE
74  if (GetVerboseLevel()>1) {
75  G4cout << "G4Decay constructor " << " Name:" << processName << G4endl;
76  }
77 #endif
78 
80 }
const G4double HighestValue
Definition: G4Decay.hh:172
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4int verboseLevel
Definition: G4Decay.hh:164
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188

Here is the call graph for this function:

G4Decay::~G4Decay ( )
virtual

Definition at line 82 of file G4Decay.cc.

83 {
84  if (pExtDecayer) {
85  delete pExtDecayer;
86  }
87 }
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181

Member Function Documentation

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

Reimplemented from G4VRestDiscreteProcess.

Reimplemented in G4DecayWithSpin.

Definition at line 191 of file G4Decay.hh.

195 {
196  return DecayIt(aTrack, aStep);
197 }
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181

Here is the call graph for this function:

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 473 of file G4Decay.cc.

477 {
478  // condition is set to "Not Forced"
479  *condition = NotForced;
480 
482  if (pTime >= 0.) {
483  fRemainderLifeTime = pTime - track.GetProperTime();
485  } else {
488  }
489  return fRemainderLifeTime;
490 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
Definition: G4Decay.cc:101
G4double GetProperTime() const
const G4DynamicParticle * GetDynamicParticle() const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
#define DBL_MIN
Definition: templates.hh:75
double G4double
Definition: G4Types.hh:76
G4double GetPreAssignedDecayProperTime() const

Here is the call graph for this function:

void G4Decay::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 176 of file G4Decay.cc.

177 {
178  return;
179 }
void G4Decay::DaughterPolarization ( const G4Track aTrack,
G4DecayProducts products 
)
protectedvirtual

Reimplemented in G4PionDecayMakeSpin.

Definition at line 383 of file G4Decay.cc.

384 {
385 }

Here is the caller graph for this function:

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

Definition at line 181 of file G4Decay.cc.

182 {
183  // The DecayIt() method returns by pointer a particle-change object.
184  // Units are expressed in GEANT4 internal units.
185 
186  // Initialize ParticleChange
187  // all members of G4VParticleChange are set to equal to
188  // corresponding member in G4Track
190 
191  // get particle
192  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
193  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
194 
195  // check if the particle is stable
196  if (aParticleDef->GetPDGStable()) return &fParticleChangeForDecay ;
197 
198 
199  //check if thePreAssignedDecayProducts exists
200  const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
201  G4bool isPreAssigned = (o_products != 0);
202  G4DecayProducts* products = 0;
203 
204  // decay table
205  G4DecayTable *decaytable = aParticleDef->GetDecayTable();
206 
207  // check if external decayer exists
208  G4bool isExtDecayer = (decaytable == 0) && (pExtDecayer !=0);
209 
210  // Error due to NO Decay Table
211  if ( (decaytable == 0) && !isExtDecayer &&!isPreAssigned ){
212  if (GetVerboseLevel()>0) {
213  G4cout << "G4Decay::DoIt : decay table not defined for ";
214  G4cout << aParticle->GetDefinition()->GetParticleName()<< G4endl;
215  }
216  G4Exception( "G4Decay::DecayIt ",
217  "DECAY101",JustWarning,
218  "Decay table is not defined");
219 
221  // Kill the parent particle
224 
226  return &fParticleChangeForDecay ;
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 = 0;
250  G4double massParent = aParticle->GetMass();
251  decaychannel = decaytable->SelectADecayChannel(massParent);
252  if ( decaychannel ==0) {
253  // decay channel not found
255  ed << "Can not determine decay channel for "
256  << aParticleDef->GetParticleName() << G4endl
257  << " mass of dynamic particle: " << massParent/GeV << " (GEV)" << G4endl
258  << " dacay table has " << decaytable->entries() << " entries" << G4endl;
259  G4double checkedmass=massParent;
260  if (massParent < 0.) {
261  checkedmass=aParticleDef->GetPDGMass();
262  ed << "Using PDG mass ("<<checkedmass/GeV << "(GeV)) in IsOKWithParentMass" << G4endl;
263  }
264  for (G4int ic =0;ic <decaytable->entries();++ic) {
265  G4VDecayChannel * dc= decaytable->GetDecayChannel(ic);
266  ed << ic << ": BR " << dc->GetBR() << ", IsOK? "
267  << dc->IsOKWithParentMass(checkedmass)
268  << ", --> ";
269  G4int ndaughters=dc->GetNumberOfDaughters();
270  for (G4int id=0;id<ndaughters;++id) {
271  if (id>0) ed << " + "; // seperator, except for first
272  ed << dc->GetDaughterName(id);
273  }
274  ed << G4endl;
275  }
276  G4Exception("G4Decay::DoIt", "DECAY003", FatalException,ed);
277  } else {
278  // execute DecayIt()
279 #ifdef G4VERBOSE
280  G4int temp = decaychannel->GetVerboseLevel();
281  if (GetVerboseLevel()>1) {
282  G4cout << "G4Decay::DoIt : selected decay channel addr:"
283  << decaychannel <<G4endl;
284  decaychannel->SetVerboseLevel(GetVerboseLevel());
285  }
286 #endif
287  products = decaychannel->DecayIt(aParticle->GetMass());
288 #ifdef G4VERBOSE
289  if (GetVerboseLevel()>1) {
290  decaychannel->SetVerboseLevel(temp);
291  }
292 #endif
293 #ifdef G4VERBOSE
294  if (GetVerboseLevel()>2) {
295  if (! products->IsChecked() ) products->DumpInfo();
296  }
297 #endif
298  }
299  }
300 
301  // get parent particle information ...................................
302  G4double ParentEnergy = aParticle->GetTotalEnergy();
303  G4double ParentMass = aParticle->GetMass();
304  if (ParentEnergy < ParentMass) {
305  if (GetVerboseLevel()>0) {
306  G4cout << "G4Decay::DoIt : Total Energy is less than its mass" << G4endl;
307  G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
308  G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
309  G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
310  G4cout << G4endl;
311  }
312  G4Exception( "G4Decay::DecayIt ",
313  "DECAY102",JustWarning,
314  "Total Energy is less than its mass");
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  {
358  // get current position of the track
359  currentPosition = aTrack.GetPosition();
360  // create a new track object
361  G4Track* secondary = new G4Track( products->PopProducts(),
362  finalGlobalTime ,
363  currentPosition );
364  // switch on good for tracking flag
365  secondary->SetGoodForTrackingFlag();
366  secondary->SetTouchableHandle(thand);
367  // add the secondary track in the List
369  }
370  delete products;
371 
372  // Kill the parent particle
375  fParticleChangeForDecay.ProposeLocalTime( finalLocalTime );
376 
377  // Clear NumberOfInteractionLengthLeft
379 
380  return &fParticleChangeForDecay ;
381 }
G4double GetBR() const
G4double GetLocalTime() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4int GetNumberOfDaughters() const
virtual G4DecayProducts * ImportDecayProducts(const G4Track &aTrack)=0
G4bool GetPDGStable() const
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VDecayChannel * GetDecayChannel(G4int index) const
virtual G4bool IsOKWithParentMass(G4double parentMass)
G4ParticleDefinition * GetDefinition() const
tuple x
Definition: test.py:50
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4int entries() const
G4DecayTable * GetDecayTable() const
virtual void Initialize(const G4Track &)
G4GLOB_DLL std::ostream G4cout
void SetVerboseLevel(G4int value)
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
void DumpInfo() const
const G4ThreeVector & GetMomentumDirection() const
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:81
static constexpr double cm
Definition: G4SIunits.hh:119
G4double GetGlobalTime() const
void AddSecondary(G4Track *aSecondary)
const G4TouchableHandle & GetTouchableHandle() const
const G4String & GetDaughterName(G4int anIndex) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int GetVerboseLevel() const
G4double GetPDGMass() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * PopProducts()
static constexpr double GeV
Definition: G4SIunits.hh:217
tuple z
Definition: test.py:28
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition: G4Decay.cc:383
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4bool IsChecked() const
G4int entries() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
#define ns
Definition: xmlparse.cc:614
void SetGoodForTrackingFlag(G4bool value=true)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Decay::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 397 of file G4Decay.cc.

398 {
399  // Clear NumberOfInteractionLengthLeft
401 
403 }
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
G4double currentInteractionLength
Definition: G4VProcess.hh:297

Here is the call graph for this function:

const G4VExtDecayer * G4Decay::GetExtDecayer ( ) const
inline

Definition at line 201 of file G4Decay.hh.

202 {
203  return pExtDecayer;
204 }
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181
G4double G4Decay::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 130 of file G4Decay.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Decay::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 101 of file G4Decay.cc.

103 {
104  // returns the mean free path in GEANT4 internal units
105  G4double meanlife;
106 
107  // get particle
108  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
109  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
110  G4double aLife = aParticleDef->GetPDGLifeTime();
111 
112  // check if the particle is stable?
113  if (aParticleDef->GetPDGStable()) {
114  //1000000 times the life time of the universe
115  meanlife = 1e24 * s;
116 
117  } else {
118  meanlife = aLife;
119  }
120 
121 #ifdef G4VERBOSE
122  if (GetVerboseLevel()>1) {
123  G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
124  }
125 #endif
126 
127  return meanlife;
128 }
const G4DynamicParticle * GetDynamicParticle() const
G4bool GetPDGStable() const
G4ParticleDefinition * GetDefinition() const
const XML_Char * s
Definition: expat.h:262
G4GLOB_DLL std::ostream G4cout
G4double GetPDGLifeTime() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define ns
Definition: xmlparse.cc:614
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Decay::GetRemainderLifeTime ( ) const
inline

Definition at line 207 of file G4Decay.hh.

208 {
209  return fRemainderLifeTime;
210 }
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4int G4Decay::GetVerboseLevel ( ) const
inline

Definition at line 188 of file G4Decay.hh.

188 { return verboseLevel; }
G4int verboseLevel
Definition: G4Decay.hh:164

Here is the caller graph for this function:

G4bool G4Decay::IsApplicable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Reimplemented in G4MuonicAtomDecay.

Definition at line 89 of file G4Decay.cc.

90 {
91  // check if the particle is stable?
92  if (aParticleType.GetPDGLifeTime() <0.0) {
93  return false;
94  } else if (aParticleType.GetPDGMass() <= 0.0*MeV) {
95  return false;
96  } else {
97  return true;
98  }
99 }
G4double GetPDGMass() const
G4double GetPDGLifeTime() const
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VRestDiscreteProcess.

Reimplemented in G4DecayWithSpin.

Definition at line 503 of file G4Decay.cc.

507 {
508  if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
509  (aTrack.GetTrackStatus() == fStopAndKill ) ){
511  return &fParticleChangeForDecay;
512  } else {
513  return DecayIt(aTrack, aStep);
514  }
515 }
G4TrackStatus GetTrackStatus() const
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181
virtual void Initialize(const G4Track &)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178

Here is the call graph for this function:

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 406 of file G4Decay.cc.

411 {
412  // condition is set to "Not Forced"
413  *condition = NotForced;
414 
415  // pre-assigned Decay time
418 
419  if (pTime < 0.) {
420  // normal case
421  if ( previousStepSize > 0.0){
422  // subtract NumberOfInteractionLengthLeft
423  SubtractNumberOfInteractionLengthLeft(previousStepSize);
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();
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 }
G4double condition(const G4ErrorSymMatrix &m)
static constexpr double perMillion
Definition: G4SIunits.hh:334
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4Decay.cc:130
G4double GetProperTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4String & GetName() const
Definition: G4Material.hh:178
void DumpInfo(G4int mode=0) const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4ParticleDefinition * GetDefinition() const
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4double GetTotalMomentum() const
G4int verboseLevel
Definition: G4Decay.hh:164
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4double GetMass() const
static constexpr double cm
Definition: G4SIunits.hh:119
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4Material * GetMaterial() const
#define DBL_MIN
Definition: templates.hh:75
G4double GetPDGLifeTime() const
#define G4endl
Definition: G4ios.hh:61
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.hh:544
double G4double
Definition: G4Types.hh:76
G4double GetPreAssignedDecayProperTime() const
#define DBL_MAX
Definition: templates.hh:83
float c_light
Definition: hepunit.py:257

Here is the call graph for this function:

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 }
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Decay::SetVerboseLevel ( G4int  value)
inline

Definition at line 185 of file G4Decay.hh.

185 { verboseLevel = value; }
G4int verboseLevel
Definition: G4Decay.hh:164
const XML_Char int const XML_Char * value
Definition: expat.h:331
void G4Decay::StartTracking ( G4Track )
virtual

Reimplemented from G4VProcess.

Definition at line 389 of file G4Decay.cc.

390 {
393 
394  fRemainderLifeTime = -1.0;
395 }
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4double currentInteractionLength
Definition: G4VProcess.hh:297

Here is the call graph for this function:

Member Data Documentation

G4ParticleChangeForDecay G4Decay::fParticleChangeForDecay
protected

Definition at line 178 of file G4Decay.hh.

G4double G4Decay::fRemainderLifeTime
protected

Definition at line 175 of file G4Decay.hh.

const G4double G4Decay::HighestValue
protected

Definition at line 172 of file G4Decay.hh.

G4VExtDecayer* G4Decay::pExtDecayer
protected

Definition at line 181 of file G4Decay.hh.

G4int G4Decay::verboseLevel
protected

Definition at line 164 of file G4Decay.hh.


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