Geant4_10
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 71228 2013-06-12 12:53:45Z 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  meanlife = DBL_MAX;
115 
116  } else {
117  meanlife = aLife;
118  }
119 
120 #ifdef G4VERBOSE
121  if (GetVerboseLevel()>1) {
122  G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
123  }
124 #endif
125 
126  return meanlife;
127 }
128 
130 {
131  // get particle
132  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
133  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
134  G4double aMass = aParticle->GetMass();
135  G4double aLife = aParticleDef->GetPDGLifeTime();
136 
137 
138  // returns the mean free path in GEANT4 internal units
139  G4double pathlength;
140  G4double aCtau = c_light * aLife;
141 
142  // check if the particle is stable?
143  if (aParticleDef->GetPDGStable()) {
144  pathlength = DBL_MAX;
145 
146  //check if the particle has very short life time ?
147  } else if (aCtau < DBL_MIN) {
148  pathlength = DBL_MIN;
149 
150  } else {
151  //calculate the mean free path
152  // by using normalized kinetic energy (= Ekin/mass)
153  G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
154  if ( rKineticEnergy > HighestValue) {
155  // gamma >> 1
156  pathlength = ( rKineticEnergy + 1.0)* aCtau;
157  } else if ( rKineticEnergy < DBL_MIN ) {
158  // too slow particle
159 #ifdef G4VERBOSE
160  if (GetVerboseLevel()>1) {
161  G4cout << "G4Decay::GetMeanFreePath() !!particle stops!!";
162  G4cout << aParticleDef->GetParticleName() << G4endl;
163  G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
164  }
165 #endif
166  pathlength = DBL_MIN;
167  } else {
168  // beta <1
169  pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
170  }
171  }
172  return pathlength;
173 }
174 
176 {
177  return;
178 }
179 
181 {
182  // The DecayIt() method returns by pointer a particle-change object.
183  // Units are expressed in GEANT4 internal units.
184 
185  // Initialize ParticleChange
186  // all members of G4VParticleChange are set to equal to
187  // corresponding member in G4Track
189 
190  // get particle
191  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
192  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
193 
194  // check if the particle is stable
195  if (aParticleDef->GetPDGStable()) return &fParticleChangeForDecay ;
196 
197 
198  //check if thePreAssignedDecayProducts exists
199  const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
200  G4bool isPreAssigned = (o_products != 0);
201  G4DecayProducts* products = 0;
202 
203  // decay table
204  G4DecayTable *decaytable = aParticleDef->GetDecayTable();
205 
206  // check if external decayer exists
207  G4bool isExtDecayer = (decaytable == 0) && (pExtDecayer !=0);
208 
209  // Error due to NO Decay Table
210  if ( (decaytable == 0) && !isExtDecayer &&!isPreAssigned ){
211  if (GetVerboseLevel()>0) {
212  G4cout << "G4Decay::DoIt : decay table not defined for ";
213  G4cout << aParticle->GetDefinition()->GetParticleName()<< G4endl;
214  }
215  G4Exception( "G4Decay::DecayIt ",
216  "DECAY101",JustWarning,
217  "Decay table is not defined");
218 
220  // Kill the parent particle
223 
225  return &fParticleChangeForDecay ;
226  }
227 
228  if (isPreAssigned) {
229  // copy decay products
230  products = new G4DecayProducts(*o_products);
231  } else if ( isExtDecayer ) {
232  // decay according to external decayer
233  products = pExtDecayer->ImportDecayProducts(aTrack);
234  } else {
235  // decay acoording to decay table
236  // choose a decay channel
237  G4VDecayChannel *decaychannel = decaytable->SelectADecayChannel();
238  if (decaychannel == 0 ){
239  // decay channel not found
240  G4Exception("G4Decay::DoIt", "DECAY003", FatalException,
241  " can not determine decay channel ");
242  } else {
243  // execute DecayIt()
244 #ifdef G4VERBOSE
245  G4int temp = decaychannel->GetVerboseLevel();
246  if (GetVerboseLevel()>1) {
247  G4cout << "G4Decay::DoIt : selected decay channel addr:" << decaychannel <<G4endl;
248  decaychannel->SetVerboseLevel(GetVerboseLevel());
249  }
250 #endif
251  products = decaychannel->DecayIt(aParticle->GetMass());
252 #ifdef G4VERBOSE
253  if (GetVerboseLevel()>1) {
254  decaychannel->SetVerboseLevel(temp);
255  }
256 #endif
257 #ifdef G4VERBOSE
258  if (GetVerboseLevel()>2) {
259  if (! products->IsChecked() ) products->DumpInfo();
260  }
261 #endif
262  }
263  }
264 
265  // get parent particle information ...................................
266  G4double ParentEnergy = aParticle->GetTotalEnergy();
267  G4double ParentMass = aParticle->GetMass();
268  if (ParentEnergy < ParentMass) {
269  if (GetVerboseLevel()>0) {
270  G4cout << "G4Decay::DoIt : Total Energy is less than its mass" << G4endl;
271  G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
272  G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
273  G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
274  G4cout << G4endl;
275  }
276  G4Exception( "G4Decay::DecayIt ",
277  "DECAY102",JustWarning,
278  "Total Energy is less than its mass");
279  ParentEnergy = ParentMass;
280  }
281 
282  G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
283 
284  //boost all decay products to laboratory frame
285  G4double energyDeposit = 0.0;
286  G4double finalGlobalTime = aTrack.GetGlobalTime();
287  G4double finalLocalTime = aTrack.GetLocalTime();
288  if (aTrack.GetTrackStatus() == fStopButAlive ){
289  // AtRest case
290  finalGlobalTime += fRemainderLifeTime;
291  finalLocalTime += fRemainderLifeTime;
292  energyDeposit += aParticle->GetKineticEnergy();
293  if (isPreAssigned) products->Boost( ParentEnergy, ParentDirection);
294  } else {
295  // PostStep case
296  if (!isExtDecayer) products->Boost( ParentEnergy, ParentDirection);
297  }
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 }
347 
349 {
350 }
351 
352 
353 
355 {
358 
359  fRemainderLifeTime = -1.0;
360 }
361 
363 {
364  // Clear NumberOfInteractionLengthLeft
366 
368 }
369 
370 
372  const G4Track& track,
373  G4double previousStepSize,
375  )
376 {
377 
378  // condition is set to "Not Forced"
379  *condition = NotForced;
380 
381  // pre-assigned Decay time
384 
385  if (pTime < 0.) {
386  // normal case
387  if ( previousStepSize > 0.0){
388  // subtract NumberOfInteractionLengthLeft
389  SubtractNumberOfInteractionLengthLeft(previousStepSize);
392  }
394  }
395  // get mean free path
396  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
397 
398 #ifdef G4VERBOSE
399  if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
400  G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
401  track.GetDynamicParticle()->DumpInfo();
402  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
403  G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
404  }
405 #endif
406 
407  G4double value;
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 }
437 
439  const G4Track& track,
441  )
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 }
456 
457 
459 {
460  pExtDecayer = val;
461 
462  // set Process Sub Type
463  if ( pExtDecayer !=0 ) {
464  SetProcessSubType(static_cast<int>(DECAY_External));
465  }
466 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
Definition: G4Decay.cc:101
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4Decay.cc:129
G4double GetLocalTime() const
virtual ~G4Decay()
Definition: G4Decay.cc:82
G4double GetProperTime() const
G4double GetKineticEnergy() const
Int_t index
Definition: macro.C:9
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
const G4double HighestValue
Definition: G4Decay.hh:172
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
virtual G4DecayProducts * ImportDecayProducts(const G4Track &aTrack)=0
const G4String & GetName() const
Definition: G4Material.hh:176
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
tuple x
Definition: test.py:50
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:354
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:180
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4double GetTotalMomentum() const
G4int verboseLevel
Definition: G4Decay.hh:164
Double_t y
Definition: plot.C:279
G4DecayTable * GetDecayTable() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4Decay.cc:371
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:438
G4double currentInteractionLength
Definition: G4VProcess.hh:297
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
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
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
G4double GetPDGMass() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * PopProducts()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
#define DBL_MIN
Definition: templates.hh:75
tuple z
Definition: test.py:28
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4Decay.cc:175
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition: G4Decay.cc:348
const XML_Char int const XML_Char * value
Definition: expat.h:331
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:458
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4double GetPreAssignedDecayProperTime() const
virtual void EndTracking()
Definition: G4Decay.cc:362
#define DBL_MAX
Definition: templates.hh:83
#define ns
Definition: xmlparse.cc:597
float c_light
Definition: hepunit.py:257
void SetGoodForTrackingFlag(G4bool value=true)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188
float perMillion
Definition: hepunit.py:241