Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4Decay Class Reference

#include <G4Decay.hh>

Inheritance diagram for G4Decay:
G4VRestDiscreteProcess G4VProcess G4DecayWithSpin G4PionDecayMakeSpin

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.

References DECAY, fParticleChangeForDecay, G4cout, G4endl, GetVerboseLevel(), G4VProcess::pParticleChange, and G4VProcess::SetProcessSubType().

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
G4Decay::~G4Decay ( )
virtual

Definition at line 82 of file G4Decay.cc.

References pExtDecayer.

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.

Definition at line 191 of file G4Decay.hh.

References DecayIt().

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 438 of file G4Decay.cc.

References DBL_MIN, fRemainderLifeTime, G4Track::GetDynamicParticle(), GetMeanLifeTime(), G4DynamicParticle::GetPreAssignedDecayProperTime(), G4Track::GetProperTime(), NotForced, and G4VProcess::theNumberOfInteractionLengthLeft.

442 {
443  // condition is set to "Not Forced"
444  *condition = NotForced;
445 
447  if (pTime >= 0.) {
448  fRemainderLifeTime = pTime - track.GetProperTime();
450  } else {
453  }
454  return fRemainderLifeTime;
455 }
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
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 348 of file G4Decay.cc.

Referenced by DecayIt().

349 {
350 }
G4VParticleChange * G4Decay::DecayIt ( const G4Track aTrack,
const G4Step aStep 
)
protectedvirtual

Reimplemented in G4DecayWithSpin.

Definition at line 181 of file G4Decay.cc.

References G4VParticleChange::AddSecondary(), G4DecayProducts::Boost(), G4VProcess::ClearNumberOfInteractionLengthLeft(), python.hepunit::cm, DaughterPolarization(), G4VDecayChannel::DecayIt(), G4DecayProducts::DumpInfo(), G4DecayProducts::entries(), FatalException, fParticleChangeForDecay, fRemainderLifeTime, fStopAndKill, fStopButAlive, G4cout, G4endl, G4Exception(), G4ParticleDefinition::GetDecayTable(), G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetLocalTime(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGStable(), G4Track::GetPosition(), G4DynamicParticle::GetPreAssignedDecayProducts(), G4DynamicParticle::GetTotalEnergy(), G4Track::GetTouchableHandle(), G4Track::GetTrackStatus(), GetVerboseLevel(), G4VDecayChannel::GetVerboseLevel(), G4VExtDecayer::ImportDecayProducts(), G4ParticleChangeForDecay::Initialize(), G4DecayProducts::IsChecked(), JustWarning, python.hepunit::MeV, ns, pExtDecayer, G4DecayProducts::PopProducts(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForDecay::ProposeLocalTime(), G4VParticleChange::ProposeTrackStatus(), G4DecayTable::SelectADecayChannel(), G4Track::SetGoodForTrackingFlag(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetTouchableHandle(), G4VDecayChannel::SetVerboseLevel(), test::x, and z.

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

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 acoording to decay table
237  // choose a decay channel
238  G4VDecayChannel *decaychannel = decaytable->SelectADecayChannel();
239  if (decaychannel == 0 ){
240  // decay channel not found
241  G4Exception("G4Decay::DoIt", "DECAY003", FatalException,
242  " can not determine decay channel ");
243  } else {
244  // execute DecayIt()
245 #ifdef G4VERBOSE
246  G4int temp = decaychannel->GetVerboseLevel();
247  if (GetVerboseLevel()>1) {
248  G4cout << "G4Decay::DoIt : selected decay channel addr:" << decaychannel <<G4endl;
249  decaychannel->SetVerboseLevel(GetVerboseLevel());
250  }
251 #endif
252  products = decaychannel->DecayIt(aParticle->GetMass());
253 #ifdef G4VERBOSE
254  if (GetVerboseLevel()>1) {
255  decaychannel->SetVerboseLevel(temp);
256  }
257 #endif
258 #ifdef G4VERBOSE
259  if (GetVerboseLevel()>2) {
260  if (! products->IsChecked() ) products->DumpInfo();
261  }
262 #endif
263  }
264  }
265 
266  // get parent particle information ...................................
267  G4double ParentEnergy = aParticle->GetTotalEnergy();
268  G4double ParentMass = aParticle->GetMass();
269  if (ParentEnergy < ParentMass) {
270  if (GetVerboseLevel()>0) {
271  G4cout << "G4Decay::DoIt : Total Energy is less than its mass" << G4endl;
272  G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
273  G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
274  G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
275  G4cout << G4endl;
276  }
277  G4Exception( "G4Decay::DecayIt ",
278  "DECAY102",JustWarning,
279  "Total Energy is less than its mass");
280  ParentEnergy = ParentMass;
281  }
282 
283  G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
284 
285  //boost all decay products to laboratory frame
286  G4double energyDeposit = 0.0;
287  G4double finalGlobalTime = aTrack.GetGlobalTime();
288  G4double finalLocalTime = aTrack.GetLocalTime();
289  if (aTrack.GetTrackStatus() == fStopButAlive ){
290  // AtRest case
291  finalGlobalTime += fRemainderLifeTime;
292  finalLocalTime += fRemainderLifeTime;
293  energyDeposit += aParticle->GetKineticEnergy();
294  if (isPreAssigned) products->Boost( ParentEnergy, ParentDirection);
295  } else {
296  // PostStep case
297  if (!isExtDecayer) products->Boost( ParentEnergy, ParentDirection);
298  }
299  // set polarization for daughter particles
300  DaughterPolarization(aTrack, products);
301 
302 
303  //add products in fParticleChangeForDecay
304  G4int numberOfSecondaries = products->entries();
305  fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
306 #ifdef G4VERBOSE
307  if (GetVerboseLevel()>1) {
308  G4cout << "G4Decay::DoIt : Decay vertex :";
309  G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
310  G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
311  G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
312  G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
313  G4cout << G4endl;
314  G4cout << "G4Decay::DoIt : decay products in Lab. Frame" << G4endl;
315  products->DumpInfo();
316  }
317 #endif
318  G4int index;
319  G4ThreeVector currentPosition;
320  const G4TouchableHandle thand = aTrack.GetTouchableHandle();
321  for (index=0; index < numberOfSecondaries; index++)
322  {
323  // get current position of the track
324  currentPosition = aTrack.GetPosition();
325  // create a new track object
326  G4Track* secondary = new G4Track( products->PopProducts(),
327  finalGlobalTime ,
328  currentPosition );
329  // switch on good for tracking flag
330  secondary->SetGoodForTrackingFlag();
331  secondary->SetTouchableHandle(thand);
332  // add the secondary track in the List
334  }
335  delete products;
336 
337  // Kill the parent particle
340  fParticleChangeForDecay.ProposeLocalTime( finalLocalTime );
341 
342  // Clear NumberOfInteractionLengthLeft
344 
345  return &fParticleChangeForDecay ;
346 }
G4double GetLocalTime() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double z
Definition: TRTMaterials.hh:39
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)
G4ParticleDefinition * GetDefinition() const
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
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
G4double GetGlobalTime() const
void AddSecondary(G4Track *aSecondary)
const G4TouchableHandle & GetTouchableHandle() const
G4VDecayChannel * SelectADecayChannel()
Definition: G4DecayTable.cc:81
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int GetVerboseLevel() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * PopProducts()
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition: G4Decay.cc:348
#define G4endl
Definition: G4ios.hh:61
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:597
void SetGoodForTrackingFlag(G4bool value=true)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188
void G4Decay::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 362 of file G4Decay.cc.

References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::currentInteractionLength.

363 {
364  // Clear NumberOfInteractionLengthLeft
366 
368 }
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
G4double currentInteractionLength
Definition: G4VProcess.hh:297
const G4VExtDecayer * G4Decay::GetExtDecayer ( ) const
inline

Definition at line 201 of file G4Decay.hh.

References pExtDecayer.

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.

References python.hepunit::c_light, DBL_MAX, DBL_MIN, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGLifeTime(), G4ParticleDefinition::GetPDGStable(), G4DynamicParticle::GetTotalMomentum(), GetVerboseLevel(), python.hepunit::GeV, and HighestValue.

Referenced by PostStepGetPhysicalInteractionLength().

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
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
G4double G4Decay::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 101 of file G4Decay.cc.

References G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4ParticleDefinition::GetPDGLifeTime(), G4ParticleDefinition::GetPDGStable(), GetVerboseLevel(), and ns.

Referenced by AtRestGetPhysicalInteractionLength().

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
const XML_Char * s
G4bool GetPDGStable() const
G4ParticleDefinition * GetDefinition() const
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:597
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188
G4double G4Decay::GetRemainderLifeTime ( ) const
inline

Definition at line 207 of file G4Decay.hh.

References fRemainderLifeTime.

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

Definition at line 188 of file G4Decay.hh.

References verboseLevel.

Referenced by DecayIt(), G4Decay(), GetMeanFreePath(), and GetMeanLifeTime().

188 { return verboseLevel; }
G4int verboseLevel
Definition: G4Decay.hh:164
G4bool G4Decay::IsApplicable ( const G4ParticleDefinition aParticleType)
virtual
G4VParticleChange * G4Decay::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 468 of file G4Decay.cc.

References DecayIt(), fParticleChangeForDecay, fStopAndKill, fStopButAlive, G4Track::GetTrackStatus(), and G4ParticleChangeForDecay::Initialize().

472 {
473  if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
474  (aTrack.GetTrackStatus() == fStopAndKill ) ){
476  return &fParticleChangeForDecay;
477  } else {
478  return DecayIt(aTrack, aStep);
479  }
480 }
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
G4double G4Decay::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 371 of file G4Decay.cc.

References python.hepunit::c_light, python.hepunit::cm, G4VProcess::currentInteractionLength, DBL_MAX, DBL_MIN, G4DynamicParticle::DumpInfo(), fRemainderLifeTime, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetMass(), G4Track::GetMaterial(), GetMeanFreePath(), G4Material::GetName(), G4ParticleDefinition::GetPDGLifeTime(), G4DynamicParticle::GetPreAssignedDecayProperTime(), G4Track::GetProperTime(), G4DynamicParticle::GetTotalMomentum(), NotForced, python.hepunit::perMillion, G4VProcess::SubtractNumberOfInteractionLengthLeft(), G4VProcess::theNumberOfInteractionLengthLeft, and verboseLevel.

376 {
377  // condition is set to "Not Forced"
378  *condition = NotForced;
379 
380  // pre-assigned Decay time
383 
384  if (pTime < 0.) {
385  // normal case
386  if ( previousStepSize > 0.0){
387  // subtract NumberOfInteractionLengthLeft
388  SubtractNumberOfInteractionLengthLeft(previousStepSize);
391  }
393  }
394  // get mean free path
395  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
396 
397 #ifdef G4VERBOSE
398  if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
399  G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
400  track.GetDynamicParticle()->DumpInfo();
401  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
402  G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
403  }
404 #endif
405 
406  G4double value;
409  //fRemainderLifeTime = theNumberOfInteractionLengthLeft*aLife;
410  } else {
411  value = DBL_MAX;
412  }
413 
414  return value;
415 
416  } else {
417  //pre-assigned Decay time case
418  // reminder proper time
419  fRemainderLifeTime = pTime - track.GetProperTime();
421 
422  G4double rvalue=0.0;
423  // use pre-assigned Decay time to determine PIL
424  if (aLife>0.0) {
425  // ordinary particle
426  rvalue = (fRemainderLifeTime/aLife)*GetMeanFreePath(track, previousStepSize, condition);
427  } else {
428  // shortlived particle
429  rvalue = c_light * fRemainderLifeTime;
430  // by using normalized kinetic energy (= Ekin/mass)
431  G4double aMass = track.GetDynamicParticle()->GetMass();
432  rvalue *= track.GetDynamicParticle()->GetTotalMomentum()/aMass;
433  }
434  return rvalue;
435  }
436 }
G4double condition(const G4ErrorSymMatrix &m)
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:176
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
G4double GetMass() const
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4Material * GetMaterial() const
#define DBL_MIN
Definition: templates.hh:75
const XML_Char int const XML_Char * value
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
float perMillion
Definition: hepunit.py:241
void G4Decay::SetExtDecayer ( G4VExtDecayer val)

Definition at line 458 of file G4Decay.cc.

References DECAY_External, pExtDecayer, and G4VProcess::SetProcessSubType().

Referenced by P6DExtDecayerPhysics::ConstructProcess().

459 {
460  pExtDecayer = val;
461 
462  // set Process Sub Type
463  if ( pExtDecayer !=0 ) {
464  SetProcessSubType(static_cast<int>(DECAY_External));
465  }
466 }
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void G4Decay::SetVerboseLevel ( G4int  value)
inline

Definition at line 185 of file G4Decay.hh.

References verboseLevel.

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

Reimplemented from G4VProcess.

Definition at line 354 of file G4Decay.cc.

References G4VProcess::currentInteractionLength, fRemainderLifeTime, and G4VProcess::ResetNumberOfInteractionLengthLeft().

355 {
358 
359  fRemainderLifeTime = -1.0;
360 }
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4double currentInteractionLength
Definition: G4VProcess.hh:297

Field Documentation

G4ParticleChangeForDecay G4Decay::fParticleChangeForDecay
protected

Definition at line 178 of file G4Decay.hh.

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

G4double G4Decay::fRemainderLifeTime
protected
const G4double G4Decay::HighestValue
protected

Definition at line 172 of file G4Decay.hh.

Referenced by GetMeanFreePath().

G4VExtDecayer* G4Decay::pExtDecayer
protected

Definition at line 181 of file G4Decay.hh.

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

G4int G4Decay::verboseLevel
protected

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