Geant4  10.02
G4Decay.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4Decay.cc 92056 2015-08-14 13:33:43Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // History: first implementation, based on object model of
34 // 2nd December 1995, G.Cosmo
35 // 7 July 1996 H.Kurashige
36 // ------------------------------------------------------------
37 // remove BuildPhysicsTable() 28 Nov. 1997 H.Kurashige
38 // change DBL_EPSIRON to DBL_MIN 14 Dec. 1997 H.Kurashige
39 // modified for new ParticleChange 12 Mar. 1998 H.Kurashige
40 // modified for "GoodForTrackingFlag" 19 June 1998 H.Kurashige
41 // rename thePhysicsTable to aPhyscisTable 2 Aug. 1998 H.Kurashige
42 // modified IsApplicable in order to protect the decay from registered
43 // to resonances 12 Dec. 1998 H.Kurashige
44 // remove G4ParticleMomentum 6 Feb. 99 H.Kurashige
45 // modified IsApplicable to activate G4Decay for resonances 1 Mar. 00 H.Kurashige
46 // Add External Decayer 23 Feb. 2001 H.Kurashige
47 // change LowestBinValue,HighestBinValue and TotBin(200) 9 Feb. 2002
48 //
49 
50 #include "G4Decay.hh"
51 
52 #include "G4PhysicalConstants.hh"
53 #include "G4SystemOfUnits.hh"
54 #include "G4DynamicParticle.hh"
55 #include "G4DecayProducts.hh"
56 #include "G4DecayTable.hh"
57 #include "G4VDecayChannel.hh"
58 #include "G4PhysicsLogVector.hh"
60 #include "G4VExtDecayer.hh"
61 
62 // constructor
63 G4Decay::G4Decay(const G4String& processName)
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 }
81 
83 {
84  if (pExtDecayer) {
85  delete pExtDecayer;
86  }
87 }
88 
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 }
100 
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 }
129 
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 }
175 
177 {
178  return;
179 }
180 
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
254  G4Exception("G4Decay::DoIt", "DECAY003", FatalException,
255  " can not determine decay channel ");
256  } else {
257  // execute DecayIt()
258 #ifdef G4VERBOSE
259  G4int temp = decaychannel->GetVerboseLevel();
260  if (GetVerboseLevel()>1) {
261  G4cout << "G4Decay::DoIt : selected decay channel addr:"
262  << decaychannel <<G4endl;
263  decaychannel->SetVerboseLevel(GetVerboseLevel());
264  }
265 #endif
266  products = decaychannel->DecayIt(aParticle->GetMass());
267 #ifdef G4VERBOSE
268  if (GetVerboseLevel()>1) {
269  decaychannel->SetVerboseLevel(temp);
270  }
271 #endif
272 #ifdef G4VERBOSE
273  if (GetVerboseLevel()>2) {
274  if (! products->IsChecked() ) products->DumpInfo();
275  }
276 #endif
277  }
278  }
279 
280  // get parent particle information ...................................
281  G4double ParentEnergy = aParticle->GetTotalEnergy();
282  G4double ParentMass = aParticle->GetMass();
283  if (ParentEnergy < ParentMass) {
284  if (GetVerboseLevel()>0) {
285  G4cout << "G4Decay::DoIt : Total Energy is less than its mass" << G4endl;
286  G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
287  G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
288  G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
289  G4cout << G4endl;
290  }
291  G4Exception( "G4Decay::DecayIt ",
292  "DECAY102",JustWarning,
293  "Total Energy is less than its mass");
294  ParentEnergy = ParentMass;
295  }
296 
297  G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
298 
299  //boost all decay products to laboratory frame
300  G4double energyDeposit = 0.0;
301  G4double finalGlobalTime = aTrack.GetGlobalTime();
302  G4double finalLocalTime = aTrack.GetLocalTime();
303  if (aTrack.GetTrackStatus() == fStopButAlive ){
304  // AtRest case
305  finalGlobalTime += fRemainderLifeTime;
306  finalLocalTime += fRemainderLifeTime;
307  energyDeposit += aParticle->GetKineticEnergy();
308  if (isPreAssigned) products->Boost( ParentEnergy, ParentDirection);
309  } else {
310  // PostStep case
311  if (!isExtDecayer) products->Boost( ParentEnergy, ParentDirection);
312  }
313  // set polarization for daughter particles
314  DaughterPolarization(aTrack, products);
315 
316 
317  //add products in fParticleChangeForDecay
318  G4int numberOfSecondaries = products->entries();
319  fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
320 #ifdef G4VERBOSE
321  if (GetVerboseLevel()>1) {
322  G4cout << "G4Decay::DoIt : Decay vertex :";
323  G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
324  G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
325  G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
326  G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
327  G4cout << G4endl;
328  G4cout << "G4Decay::DoIt : decay products in Lab. Frame" << G4endl;
329  products->DumpInfo();
330  }
331 #endif
332  G4int index;
333  G4ThreeVector currentPosition;
334  const G4TouchableHandle thand = aTrack.GetTouchableHandle();
335  for (index=0; index < numberOfSecondaries; index++)
336  {
337  // get current position of the track
338  currentPosition = aTrack.GetPosition();
339  // create a new track object
340  G4Track* secondary = new G4Track( products->PopProducts(),
341  finalGlobalTime ,
342  currentPosition );
343  // switch on good for tracking flag
344  secondary->SetGoodForTrackingFlag();
345  secondary->SetTouchableHandle(thand);
346  // add the secondary track in the List
348  }
349  delete products;
350 
351  // Kill the parent particle
354  fParticleChangeForDecay.ProposeLocalTime( finalLocalTime );
355 
356  // Clear NumberOfInteractionLengthLeft
358 
359  return &fParticleChangeForDecay ;
360 }
361 
363 {
364 }
365 
366 
367 
369 {
372 
373  fRemainderLifeTime = -1.0;
374 }
375 
377 {
378  // Clear NumberOfInteractionLengthLeft
380 
382 }
383 
384 
386  const G4Track& track,
387  G4double previousStepSize,
389  )
390 {
391  // condition is set to "Not Forced"
392  *condition = NotForced;
393 
394  // pre-assigned Decay time
397 
398  if (pTime < 0.) {
399  // normal case
400  if ( previousStepSize > 0.0){
401  // subtract NumberOfInteractionLengthLeft
402  SubtractNumberOfInteractionLengthLeft(previousStepSize);
405  }
407  }
408  // get mean free path
409  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
410 
411 #ifdef G4VERBOSE
412  if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
413  G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
414  track.GetDynamicParticle()->DumpInfo();
415  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
416  G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
417  }
418 #endif
419 
420  G4double value;
423  //fRemainderLifeTime = theNumberOfInteractionLengthLeft*aLife;
424  } else {
425  value = DBL_MAX;
426  }
427 
428  return value;
429 
430  } else {
431  //pre-assigned Decay time case
432  // reminder proper time
433  fRemainderLifeTime = pTime - track.GetProperTime();
435 
436  G4double rvalue=0.0;
437  // use pre-assigned Decay time to determine PIL
438  if (aLife>0.0) {
439  // ordinary particle
440  rvalue = (fRemainderLifeTime/aLife)*GetMeanFreePath(track, previousStepSize, condition);
441  } else {
442  // shortlived particle
443  rvalue = c_light * fRemainderLifeTime;
444  // by using normalized kinetic energy (= Ekin/mass)
445  G4double aMass = track.GetDynamicParticle()->GetMass();
446  rvalue *= track.GetDynamicParticle()->GetTotalMomentum()/aMass;
447  }
448  return rvalue;
449  }
450 }
451 
453  const G4Track& track,
455  )
456 {
457  // condition is set to "Not Forced"
458  *condition = NotForced;
459 
461  if (pTime >= 0.) {
462  fRemainderLifeTime = pTime - track.GetProperTime();
464  } else {
467  }
468  return fRemainderLifeTime;
469 }
470 
471 
473 {
474  pExtDecayer = val;
475 
476  // set Process Sub Type
477  if ( pExtDecayer !=0 ) {
478  SetProcessSubType(static_cast<int>(DECAY_External));
479  }
480 }
481 
483  const G4Track& aTrack,
484  const G4Step& aStep
485  )
486 {
487  if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
488  (aTrack.GetTrackStatus() == fStopAndKill ) ){
490  return &fParticleChangeForDecay;
491  } else {
492  return DecayIt(aTrack, aStep);
493  }
494 }
495 
static const double cm
Definition: G4SIunits.hh:118
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
Definition: G4Decay.cc:101
static const double MeV
Definition: G4SIunits.hh:211
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4Decay.cc:130
G4double GetLocalTime() const
virtual ~G4Decay()
Definition: G4Decay.cc:82
G4double GetProperTime() const
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
const G4double HighestValue
Definition: G4Decay.hh:172
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
G4double z
Definition: TRTMaterials.hh:39
virtual G4DecayProducts * ImportDecayProducts(const G4Track &aTrack)=0
G4bool GetPDGStable() const
const G4String & GetName() const
Definition: G4Material.hh:178
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
void DumpInfo(G4int mode=0) const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4ParticleDefinition * GetDefinition() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:482
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:181
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
virtual void StartTracking(G4Track *)
Definition: G4Decay.cc:368
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4double GetTotalMomentum() const
G4int verboseLevel
Definition: G4Decay.hh:164
static const double s
Definition: G4SIunits.hh:168
G4DecayTable * GetDecayTable() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4Decay.cc:385
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
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
Definition: G4Decay.cc:452
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:81
G4double currentInteractionLength
Definition: G4VProcess.hh:297
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double GeV
Definition: G4SIunits.hh:214
Definition: G4Step.hh:76
G4Decay(const G4String &processName="Decay")
Definition: G4Decay.cc:63
G4double GetGlobalTime() const
void AddSecondary(G4Track *aSecondary)
const G4TouchableHandle & GetTouchableHandle() const
G4Material * GetMaterial() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double perMillion
Definition: G4SIunits.hh:331
G4int GetVerboseLevel() const
G4double GetPDGMass() const
void SetNumberOfSecondaries(G4int totSecondaries)
const G4double x[NPOINTSGL]
G4DynamicParticle * PopProducts()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
#define DBL_MIN
Definition: templates.hh:75
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4Decay.cc:176
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition: G4Decay.cc:362
G4double GetPDGLifeTime() const
#define G4endl
Definition: G4ios.hh:61
G4bool IsChecked() const
G4int entries() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.hh:544
void SetExtDecayer(G4VExtDecayer *)
Definition: G4Decay.cc:472
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4double GetPreAssignedDecayProperTime() const
virtual void EndTracking()
Definition: G4Decay.cc:376
#define DBL_MAX
Definition: templates.hh:83
#define ns
Definition: xmlparse.cc:614
void SetGoodForTrackingFlag(G4bool value=true)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188