Geant4_10
G4RadioactiveDecay.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 // MODULE: G4RadioactiveDecay.cc
27 //
28 // Author: F Lei & P R Truscott
29 // Organisation: DERA UK
30 // Customer: ESA/ESTEC, NOORDWIJK
31 // Contract: 12115/96/JG/NL Work Order No. 3
32 //
33 // Documentation avaialable at http://www.space.dera.gov.uk/space_env/rdm.html
34 // These include:
35 // User Requirement Document (URD)
36 // Software Specification Documents (SSD)
37 // Software User Manual (SUM)
38 // Technical Note (TN) on the physics and algorithms
39 //
40 // The test and example programs are not included in the public release of
41 // G4 but they can be downloaded from
42 // http://www.space.qinetiq.com/space_env/rdm.html
43 //
44 // CHANGE HISTORY
45 // --------------
46 // 03 Oct 2012, V. Ivanchenko removed internal table for mean free path
47 // similar to what is done for as G4Decay
48 // 10 July 2012, L. Desorgher
49 // -In LoadDecayTable:
50 // Add LoadedNuclei.push_back(theParentNucleus.GetParticleName());
51 // also for the case where user data files are used. Correction for bug
52 // 1324. Changes proposed by Joa L.
53 //
54 //
55 // 01 May 2012, L. Desorgher
56 // -Force the reading of user file to theIsotopeTable
57 // -Merge the development by Fan Lei for activation computation
58 //
59 // 17 Oct 2011, L. Desorgher
60 // -Add possibility for the user to load its own decay file.
61 // -Set halflifethreshold negative by default to allow the tracking of all
62 // excited nuclei resulting from a radioactive decay
63 //
64 // 01 June 2011, M. Kelsey -- Add directional biasing interface to allow for
65 // "collimation" of decay daughters.
66 // 16 February 2006, V.Ivanchenko fix problem in IsApplicable connected with
67 // 8.0 particle design
68 // 18 October 2002, F. Lei
69 // in the case of beta decay, added a check of the end-energy
70 // to ensure it is > 0.
71 // ENSDF occationally have beta decay entries with zero energies
72 //
73 // 27 Sepetember 2001, F. Lei
74 // verboselevel(0) used in constructor
75 //
76 // 01 November 2000, F.Lei
77 // added " ee = e0 +1. ;" as line 763
78 // tagged as "radiative_decay-V02-00-02"
79 // 28 October 2000, F Lei
80 // added fast beta decay mode. Many files have been changed.
81 // tagged as "radiative_decay-V02-00-01"
82 //
83 // 25 October 2000, F Lei, DERA UK
84 // 1) line 1185 added 'const' to work with tag "Track-V02-00-00"
85 // tagged as "radiative_decay-V02-00-00"
86 // 14 April 2000, F Lei, DERA UK
87 // 0.b.4 release. Changes are:
88 // 1) Use PhotonEvaporation instead of DiscreteGammaDeexcitation
89 // 2) VR: Significant efficiency improvement
90 //
91 // 29 February 2000, P R Truscott, DERA UK
92 // 0.b.3 release.
93 //
95 //
96 #include "G4RadioactiveDecay.hh"
98 
99 #include "G4SystemOfUnits.hh"
100 #include "G4DynamicParticle.hh"
101 #include "G4DecayProducts.hh"
102 #include "G4DecayTable.hh"
104 #include "G4ITDecayChannel.hh"
106 #include "G4BetaPlusDecayChannel.hh"
107 #include "G4KshellECDecayChannel.hh"
108 #include "G4LshellECDecayChannel.hh"
109 #include "G4MshellECDecayChannel.hh"
110 #include "G4AlphaDecayChannel.hh"
111 #include "G4VDecayChannel.hh"
112 #include "G4RadioactiveDecayMode.hh"
113 #include "G4Ions.hh"
114 #include "G4IonTable.hh"
115 #include "G4RIsotopeTable.hh"
116 #include "G4BetaDecayType.hh"
117 #include "G4BetaDecayCorrections.hh"
118 #include "Randomize.hh"
119 #include "G4LogicalVolumeStore.hh"
120 #include "G4NuclearLevelManager.hh"
121 #include "G4NuclearLevelStore.hh"
122 #include "G4ThreeVector.hh"
123 #include "G4Electron.hh"
124 #include "G4Positron.hh"
125 #include "G4Neutron.hh"
126 #include "G4Gamma.hh"
127 #include "G4Alpha.hh"
128 
129 #include "G4HadronicProcessType.hh"
130 #include "G4LossTableManager.hh"
131 #include "G4VAtomDeexcitation.hh"
132 
133 #include <vector>
134 #include <sstream>
135 #include <algorithm>
136 #include <fstream>
137 
138 using namespace CLHEP;
139 
140 const G4double G4RadioactiveDecay::levelTolerance = 2.0*keV;
141 const G4ThreeVector G4RadioactiveDecay::origin(0.,0.,0.);
142 #ifdef G4MULTITHREADED
143 DecayTableMap* G4RadioactiveDecay::master_dkmap = 0;
144 #endif
145 
147  : G4VRestDiscreteProcess(processName, fDecay), isInitialised(false),
148  forceDecayDirection(0.,0.,0.), forceDecayHalfAngle(0.*deg), verboseLevel(0)
149 {
150 #ifdef G4VERBOSE
151  if (GetVerboseLevel() > 1) {
152  G4cout << "G4RadioactiveDecay constructor: processName = " << processName
153  << G4endl;
154  }
155 #endif
156 
158 
159  theRadioactiveDecaymessenger = new G4RadioactiveDecaymessenger(this);
160  pParticleChange = &fParticleChangeForRadDecay;
161  theIsotopeTable = new G4RIsotopeTable();
162 
163  // Regsiter the isotope table to the ion table.
164  // Although we are touching the ion table, which is shared, we are not
165  // adding particles in this operation. We can therefore do the registration
166  // for each instance of the RDM process and do not need to restrict it to the
167  // master process. It's possible that future, more optimized versions of
168  // G4IonTable will require to test for master.
169 
171  theIonTable->RegisterIsotopeTable(theIsotopeTable);
172 
173  // Reset the list of user defined data files
174  theUserRadioactiveDataFiles.clear();
175 
176  // Instantiate the map of decay tables
177 #ifdef G4MULTITHREADED
178  if(!master_dkmap) master_dkmap = new DecayTableMap;
179 #endif
180  dkmap = new DecayTableMap;
181 
182  // Apply default values.
183  NSourceBin = 1;
184  SBin[0] = 0.* s;
185  SBin[1] = 1.* s;
186  SProfile[0] = 1.;
187  SProfile[1] = 0.;
188  NDecayBin = 1;
189  DBin[0] = 0. * s ;
190  DBin[1] = 1. * s;
191  DProfile[0] = 1.;
192  DProfile[1] = 0.;
193  decayWindows[0] = 0;
195  theRadioactivityTables.push_back(rTable);
196  NSplit = 1;
197  AnalogueMC = true ;
198  FBeta = false ;
199  BRBias = true ;
200  applyICM = true ;
201  applyARM = true ;
202  halflifethreshold = nanosecond;
203 
204  // RDM applies to xall logical volumes as default
205  isAllVolumesMode = true;
207 }
208 
209 
211 {
212  delete theRadioactiveDecaymessenger;
213  delete theIsotopeTable;
214 }
215 
216 
218 {
219  // All particles other than G4Ions, are rejected by default
220  if (((const G4Ions*)(&aParticle))->GetExcitationEnergy() > 0.) {return true;}
221  if (aParticle.GetParticleName() == "GenericIon") {
222  return true;
223  } else if (!(aParticle.GetParticleType() == "nucleus")
224  || aParticle.GetPDGLifeTime() < 0. ) {
225  return false;
226  }
227 
228  // Determine whether the nuclide falls into the correct A and Z range
229  G4int A = ((const G4Ions*) (&aParticle))->GetAtomicMass();
230  G4int Z = ((const G4Ions*) (&aParticle))->GetAtomicNumber();
231  if (A > theNucleusLimits.GetAMax() || A < theNucleusLimits.GetAMin())
232  {return false;}
233  else if (Z > theNucleusLimits.GetZMax() || Z < theNucleusLimits.GetZMin())
234  {return false;}
235  return true;
236 }
237 
239 {
240  G4String key = aNucleus->GetParticleName();
241  DecayTableMap::iterator table_ptr = dkmap->find(key);
242 
243  G4DecayTable* theDecayTable = 0;
244  if (table_ptr == dkmap->end() ) { // If table not there,
245  theDecayTable = LoadDecayTable(*aNucleus); // load from file and
246  (*dkmap)[key] = theDecayTable; // store in library
247  } else {
248  theDecayTable = table_ptr->second;
249  }
250 
251  return theDecayTable;
252 }
253 
254 
256 {
257  G4LogicalVolumeStore* theLogicalVolumes;
258  G4LogicalVolume* volume;
259  theLogicalVolumes = G4LogicalVolumeStore::GetInstance();
260  for (size_t i = 0; i < theLogicalVolumes->size(); i++) {
261  volume=(*theLogicalVolumes)[i];
262  if (volume->GetName() == aVolume) {
263  ValidVolumes.push_back(aVolume);
264  std::sort(ValidVolumes.begin(), ValidVolumes.end());
265  // sort need for performing binary_search
266 #ifdef G4VERBOSE
267  if (GetVerboseLevel()>0)
268  G4cout << " RDM Applies to : " << aVolume << G4endl;
269 #endif
270  } else if(i == theLogicalVolumes->size()) {
271  G4cerr << "SelectAVolume: "<< aVolume
272  << " is not a valid logical volume name" << G4endl;
273  }
274  }
275 }
276 
277 
279 {
280  G4LogicalVolumeStore* theLogicalVolumes;
281  G4LogicalVolume* volume;
282  theLogicalVolumes = G4LogicalVolumeStore::GetInstance();
283  for (size_t i = 0; i < theLogicalVolumes->size(); i++){
284  volume=(*theLogicalVolumes)[i];
285  if (volume->GetName() == aVolume) {
286  std::vector<G4String>::iterator location;
287  location = std::find(ValidVolumes.begin(),ValidVolumes.end(),aVolume);
288  if (location != ValidVolumes.end()) {
289  ValidVolumes.erase(location);
290  std::sort(ValidVolumes.begin(), ValidVolumes.end());
291  isAllVolumesMode =false;
292  } else {
293  G4cerr << " DeselectVolume:" << aVolume << " is not in the list "
294  << G4endl;
295  }
296 #ifdef G4VERBOSE
297  if (GetVerboseLevel() > 0)
298  G4cout << " DeselectVolume: " << aVolume << " is removed from list "
299  << G4endl;
300 #endif
301  } else if (i == theLogicalVolumes->size()) {
302  G4cerr << " DeselectVolume:" << aVolume
303  << "is not a valid logical volume name" << G4endl;
304  }
305  }
306 }
307 
308 
310 {
311  G4LogicalVolumeStore* theLogicalVolumes;
312  G4LogicalVolume* volume;
313  theLogicalVolumes = G4LogicalVolumeStore::GetInstance();
314  ValidVolumes.clear();
315 #ifdef G4VERBOSE
316  if (GetVerboseLevel()>0)
317  G4cout << " RDM Applies to all Volumes" << G4endl;
318 #endif
319  for (size_t i = 0; i < theLogicalVolumes->size(); i++){
320  volume = (*theLogicalVolumes)[i];
321  ValidVolumes.push_back(volume->GetName());
322 #ifdef G4VERBOSE
323  if (GetVerboseLevel()>0)
324  G4cout << " RDM Applies to Volume " << volume->GetName() << G4endl;
325 #endif
326  }
327  std::sort(ValidVolumes.begin(), ValidVolumes.end());
328  // sort needed in order to allow binary_search
329  isAllVolumesMode=true;
330 }
331 
332 
334 {
335  ValidVolumes.clear();
336  isAllVolumesMode=false;
337 #ifdef G4VERBOSE
338  if (GetVerboseLevel() > 0) G4cout << "RDM removed from all volumes" << G4endl;
339 #endif
340 }
341 
342 
343 G4bool
345 {
346  // Check whether the radioactive decay rates table for the ion has already
347  // been calculated.
348  G4String aParticleName = aParticle.GetParticleName();
349  for (size_t i = 0; i < theDecayRateTableVector.size(); i++) {
350  if (theDecayRateTableVector[i].GetIonName() == aParticleName) return true;
351  }
352  return false;
353 }
354 
355 // GetDecayRateTable
356 // retrieve the decayratetable for the specified aParticle
357 
358 void
360 {
361  G4String aParticleName = aParticle.GetParticleName();
362 
363  for (size_t i = 0; i < theDecayRateTableVector.size(); i++) {
364  if (theDecayRateTableVector[i].GetIonName() == aParticleName) {
365  theDecayRateVector = theDecayRateTableVector[i].GetItsRates();
366  }
367  }
368 #ifdef G4VERBOSE
369  if (GetVerboseLevel() > 0) {
370  G4cout << "The DecayRate Table for " << aParticleName << " is selected."
371  << G4endl;
372  }
373 #endif
374 }
375 
376 // GetTaoTime performs the convolution of the source time profile function
377 // with the decay constants in the decay chain.
379 {
380  long double taotime = 0.L;
381  G4int nbin;
382  if ( t > SBin[NSourceBin]) {
383  nbin = NSourceBin;}
384  else {
385  nbin = 0;
386  while (t > SBin[nbin]) nbin++;
387  nbin--;}
388  long double lt = t ;
389  long double ltao = tao;
390 
391  if (nbin > 0) {
392  for (G4int i = 0; i < nbin; i++) {
393  taotime += (long double)SProfile[i] *
394  (std::exp(-(lt-(long double)SBin[i+1])/ltao)-std::exp(-(lt-(long double)SBin[i])/ltao));
395  }
396  }
397  taotime += (long double)SProfile[nbin] * (1.L-std::exp(-(lt-(long double)SBin[nbin])/ltao));
398  if (taotime < 0.) {
399  G4cout <<" Tao time =: " <<taotime << " reset to zero!"<<G4endl;
400  G4cout <<" t = " << t <<" tao = " <<tao <<G4endl;
401  G4cout << SBin[nbin] << " " <<SBin[0] << G4endl;
402  taotime = 0.;
403  }
404 #ifdef G4VERBOSE
405  if (GetVerboseLevel()>1)
406  {G4cout <<" Tao time: " <<taotime <<G4endl;}
407 #endif
408  return (G4double)taotime ;
409 }
410 
411 /*
412 // Other implementation tests to avoid use of long double
413 G4double G4RadioactiveDecay::GetTaoTime(const G4double t, const G4double tao)
414 {
415  long double taotime =0.L;
416  G4int nbin;
417  if ( t > SBin[NSourceBin]) {
418  nbin = NSourceBin;}
419  else {
420  nbin = 0;
421  while (t > SBin[nbin]) nbin++;
422  nbin--;}
423  long double lt = t ;
424  long double ltao = tao;
425  long double factor,factor1,dt1,dt;
426  if (nbin > 0) {
427  for (G4int i = 0; i < nbin; i++)
428  { long double s1=SBin[i];
429  long double s2=SBin[i+1];
430  dt1=(s2-s1)/ltao;
431  if (dt1 <50.) {
432  factor1=std::exp(dt1)-1.;
433  if (factor1<dt1) factor1 =dt1;
434  dt=(lt-s1)/ltao;
435  factor=std::exp(-dt);
436  }
437  else {
438  factor1=1.-std::exp(-dt1);
439  dt=(lt-s2)/ltao;
440  factor=std::exp(-dt);
441  }
442  G4cout<<(long double) SProfile[i] *factor*factor1<<'\t'<<std::endl;
443  long double test = (long double)SProfile[i] * (std::exp(-(lt-(long double)SBin[i+1])/ltao)-std::exp(-(lt-(long double)SBin[i])/ltao));
444  G4cout<<test<<std::endl;
445  taotime += (long double) SProfile[i] *factor*factor1;
446  }
447  }
448  long double s=SBin[nbin];
449  dt1=(lt-s)/ltao;
450  factor=1.-std::exp(-dt1);
451  taotime += (long double) SProfile[nbin] *factor;
452  if (taotime < 0.) {
453  G4cout <<" Tao time =: " <<taotime << " reset to zero!"<<G4endl;
454  G4cout <<" t = " << t <<" tao = " <<tao <<G4endl;
455  G4cout << SBin[nbin] << " " <<SBin[0] << G4endl;
456  taotime = 0.;
457  }
458 #ifdef G4VERBOSE
459  if (GetVerboseLevel()>1)
460  {G4cout <<" Tao time: " <<taotime <<G4endl;}
461 #endif
462  return (G4double)taotime ;
463 }
464 
465 
466 
467 G4double G4RadioactiveDecay::GetTaoTime(const G4double t, const G4double tao)
468 {
469  G4double taotime =0.;
470  G4int nbin;
471  if ( t > SBin[NSourceBin]) {
472  nbin = NSourceBin;}
473  else {
474  nbin = 0;
475  while (t > SBin[nbin]) nbin++;
476  nbin--;}
477  G4double lt = t ;
478  G4double ltao = tao;
479  G4double factor,factor1,dt1,dt;
480 
481  if (nbin > 0) {
482  for (G4int i = 0; i < nbin; i++)
483  { dt1=(SBin[i+1]-SBin[i])/ltao;
484  if (dt1 <50.) {
485  factor1=std::exp(dt1)-1.;
486  if (factor1<dt1) factor1 =dt1;
487  dt=(lt-SBin[i])/ltao;
488  factor=std::exp(-(lt-SBin[i])/ltao);
489  G4cout<<factor<<'\t'<<factor1<<std::endl;
490  }
491  else {
492  factor1=1.-std::exp(-dt1);
493  factor=std::exp(-(lt-SBin[i+1])/ltao);
494  }
495  G4cout<<factor<<'\t'<<factor1<<std::endl;
496  taotime += SProfile[i] *factor*factor1;
497  G4cout<<taotime<<std::endl;
498  }
499  }
500  dt1=(lt-SBin[nbin])/ltao;
501  factor=1.-std::exp(-dt1);
502  if (factor<(dt1-0.5*dt1*dt1)) factor =dt1-0.5*dt1*dt1;
503 
504 
505  taotime += SProfile[nbin] *factor;
506  G4cout<<factor<<'\t'<<taotime<<std::endl;
507  if (taotime < 0.) {
508  G4cout <<" Tao time =: " <<taotime << " reset to zero!"<<G4endl;
509  G4cout <<" t = " << t <<" tao = " <<tao <<G4endl;
510  G4cout << SBin[nbin] << " " <<SBin[0] << G4endl;
511  taotime = 0.;
512  }
513 
514 #ifdef G4VERBOSE
515  if (GetVerboseLevel()>1)
516  {G4cout <<" Tao time: " <<taotime <<G4endl;}
517 #endif
518  return (G4double)taotime ;
519 }
520 */
522 // //
523 // GetDecayTime //
524 // Randomly select a decay time for the decay process, following the //
525 // supplied decay time bias scheme. //
526 // //
528 
530 {
531  G4double decaytime = 0.;
532  G4double rand = G4UniformRand();
533  G4int i = 0;
534  while ( DProfile[i] < rand) i++;
535  rand = G4UniformRand();
536  decaytime = DBin[i] + rand*(DBin[i+1]-DBin[i]);
537 #ifdef G4VERBOSE
538  if (GetVerboseLevel()>1)
539  {G4cout <<" Decay time: " <<decaytime/s <<"[s]" <<G4endl;}
540 #endif
541  return decaytime;
542 }
543 
544 
546 {
547  G4int i = 0;
548  while ( aDecayTime > DBin[i] ) i++;
549  return i;
550 }
551 
553 // //
554 // GetMeanLifeTime (required by the base class) //
555 // //
557 
560 {
561  // For varience reduction implementation the time is set to 0 so as to
562  // force the particle to decay immediately.
563  // In analogueMC mode it return the particle's mean-life.
564 
565  G4double meanlife = 0.;
566  if (AnalogueMC) {
567  const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
568  G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
569  G4double theLife = theParticleDef->GetPDGLifeTime();
570 
571 #ifdef G4VERBOSE
572  if (GetVerboseLevel() > 2) {
573  G4cout << "G4RadioactiveDecay::GetMeanLifeTime() " << G4endl;
574  G4cout << "KineticEnergy: " << theParticle->GetKineticEnergy()/GeV
575  << " GeV, Mass: " << theParticle->GetMass()/GeV
576  << " GeV, Life time: " << theLife/ns << " ns " << G4endl;
577  }
578 #endif
579  if (theParticleDef->GetPDGStable()) {meanlife = DBL_MAX;}
580  else if (theLife < 0.0) {meanlife = DBL_MAX;}
581  else {meanlife = theLife;}
582  // Set meanlife to zero for excited istopes which are not in the
583  // RDM database
584  if (((const G4Ions*)(theParticleDef))->GetExcitationEnergy() > 0. &&
585  meanlife == DBL_MAX) {meanlife = 0.;}
586  }
587 #ifdef G4VERBOSE
588  if (GetVerboseLevel() > 1)
589  G4cout << " mean life time: " << meanlife/s << " s " << G4endl;
590 #endif
591 
592  return meanlife;
593 }
594 
596 // //
597 // GetMeanFreePath for decay in flight //
598 // //
600 
603 {
604  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
605  G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
606  G4double tau = aParticleDef->GetPDGLifeTime();
607  G4double aMass = aParticle->GetMass();
608 
609 #ifdef G4VERBOSE
610  if (GetVerboseLevel() > 2) {
611  G4cout << "G4RadioactiveDecay::GetMeanFreePath() " << G4endl;
612  G4cout << " KineticEnergy: " << aParticle->GetKineticEnergy()/GeV
613  << " GeV, Mass: " << aMass/GeV << " GeV, tau: " << tau << " ns "
614  << G4endl;
615  }
616 #endif
617  G4double pathlength = DBL_MAX;
618  if (tau != -1) {
619  // Ion can decay
620 
621  if (tau < 0.0) {
623  ed << "Ion has negative lifetime " << tau
624  << " but is not stable. Setting mean free path to DBL_MAX" << G4endl;
625  G4Exception("G4RadioactiveDecay::GetMeanFreePath()", "HAD_RDM_011",
626  JustWarning, ed);
627  pathlength = DBL_MAX;
628 
629  } else {
630  // Calculate mean free path
631  G4double betaGamma = aParticle->GetTotalMomentum()/aMass;
632  pathlength = c_light*tau*betaGamma;
633 
634  if (pathlength < DBL_MIN) {
635  pathlength = DBL_MIN;
636 #ifdef G4VERBOSE
637  if (GetVerboseLevel() > 2) {
638  G4cout << "G4Decay::GetMeanFreePath: "
639  << aParticleDef->GetParticleName()
640  << " stops, kinetic energy = "
641  << aParticle->GetKineticEnergy()/keV <<" keV " << G4endl;
642  }
643 #endif
644  }
645  }
646  }
647 
648 #ifdef G4VERBOSE
649  if (GetVerboseLevel() > 1) {
650  G4cout << "mean free path: "<< pathlength/m << " m" << G4endl;
651  }
652 #endif
653  return pathlength;
654 }
655 
657 // //
658 // BuildPhysicsTable - initialisation of atomic de-excitation //
659 // //
661 
663 {
664  if (!isInitialised) {
665  isInitialised = true;
667  G4VAtomDeexcitation* p = theManager->AtomDeexcitation();
668  if(p) { p->InitialiseAtomicDeexcitation(); }
669  }
670 }
671 
673 // //
674 // LoadDecayTable loads the decay scheme from the RadioactiveDecay database //
675 // for the parent nucleus. //
676 // //
678 
679 #ifdef G4MULTITHREADED
680 #include "G4AutoLock.hh"
681 G4Mutex G4RadioactiveDecay::radioactiveDecayMutex = G4MUTEX_INITIALIZER;
682 #endif
683 
686 {
687  // Generate input data file name using Z and A of the parent nucleus
688  // file containing radioactive decay data.
689  G4int A = ((const G4Ions*)(&theParentNucleus))->GetAtomicMass();
690  G4int Z = ((const G4Ions*)(&theParentNucleus))->GetAtomicNumber();
691  G4int lvl = ((const G4Ions*)(&theParentNucleus))->GetIsomerLevel();
692  G4DecayTable* theDecayTable = 0;
693 
694 #ifdef G4MULTITHREADED
695  G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
696 
697  G4String key = theParentNucleus.GetParticleName();
698  DecayTableMap::iterator master_table_ptr = master_dkmap->find(key);
699 
700  if (master_table_ptr != master_dkmap->end() ) { // If table is there
701  return master_table_ptr->second;
702  }
703 #endif
704 
705  // Create and initialise variables used in the method.
706  theDecayTable = new G4DecayTable();
707 
708  //Check if data have been provided by the user
709  G4String file= theUserRadioactiveDataFiles[1000*A+Z];
710 
711  if (file =="") {
712  if (!getenv("G4RADIOACTIVEDATA") ) {
713  G4cout << "Please setenv G4RADIOACTIVEDATA to point to the radioactive decay data files."
714  << G4endl;
715  throw G4HadronicException(__FILE__, __LINE__, " Please setenv G4RADIOACTIVEDATA to point to the radioactive decay data files.");
716  }
717  G4String dirName = getenv("G4RADIOACTIVEDATA");
718 
719  std::ostringstream os;
720  os <<dirName <<"/z" <<Z <<".a" <<A <<'\0';
721  file = os.str();
722  }
723 
724  std::ifstream DecaySchemeFile(file);
725 
726  G4bool found(false);
727  if (DecaySchemeFile) {
728  // Initialise variables used for reading in radioactive decay data.
729  G4int nMode = 7;
730  G4bool modeFirstRecord[7];
731  G4double modeTotalBR[7] = {0.0};
732  G4double modeSumBR[7];
733  for (G4int i = 0; i < nMode; i++) {
734  modeFirstRecord[i] = true;
735  modeSumBR[i] = 0.0;
736  }
737 
738  G4bool complete(false);
739  char inputChars[100]={' '};
740  G4String inputLine;
741  G4String recordType("");
742  G4RadioactiveDecayMode theDecayMode;
743  G4double a(0.0);
744  G4double b(0.0);
745  G4double c(0.0);
746  G4int levelCounter = 0;
747  G4BetaDecayType betaType(allowed);
748  G4double e0;
749 
750  // Loop through each data file record until you identify the decay
751  // data relating to the nuclide of concern.
752 
753  while (!complete && !DecaySchemeFile.getline(inputChars, 100).eof()) {
754  inputLine = inputChars;
755  inputLine = inputLine.strip(1);
756  if (inputChars[0] != '#' && inputLine.length() != 0) {
757  std::istringstream tmpStream(inputLine);
758 
759  if (inputChars[0] == 'P') {
760  // Nucleus is a parent type. Check excitation level to see if it
761  // matches that of theParentNucleus
762  tmpStream >> recordType >> a >> b;
763  if (found) {
764  complete = true;
765 // else {found = (std::abs(a*keV - E) < levelTolerance);}
766  } else {
767  found = (levelCounter == lvl);
768  }
769  levelCounter++;
770 
771  } else if (found) {
772  // The right part of the radioactive decay data file has been found. Search
773  // through it to determine the mode of decay of the subsequent records.
774  if (inputChars[0] == 'W') {
775 #ifdef G4VERBOSE
776  if (GetVerboseLevel() > 0) {
777  // a comment line identified and print out the message
778  G4cout << " Warning in G4RadioactiveDecay::LoadDecayTable " << G4endl;
779  G4cout << " In data file " << file << G4endl;
780  G4cout << " " << inputLine << G4endl;
781  }
782 #endif
783  } else {
784  tmpStream >> theDecayMode >> a >> b >> c >> betaType;
785 
786  // Allowed transitions are the default. Forbidden transitions are
787  // indicated in the last column.
788  if (inputLine.length() < 80) betaType = allowed;
789  a /= 1000.;
790  c /= 1000.;
791 
792  switch (theDecayMode) {
793 
794  case IT: // Isomeric transition
795  {
796  G4ITDecayChannel* anITChannel =
798  (const G4Ions*)& theParentNucleus, b);
799  anITChannel->SetICM(applyICM);
800  anITChannel->SetARM(applyARM);
801  anITChannel->SetHLThreshold(halflifethreshold);
802  theDecayTable->Insert(anITChannel);
803  }
804  break;
805 
806  case BetaMinus:
807  {
808  if (modeFirstRecord[1]) {
809  modeFirstRecord[1] = false;
810  modeTotalBR[1] = b;
811  } else {
812  if (c > 0.) {
813  e0 = c*MeV/0.511;
814  G4BetaDecayCorrections corrections(Z+1, A);
815 
816  // array to store function shape
817  G4int npti = 100;
818  G4double* pdf = new G4double[npti];
819 
820  G4double e; // Total electron energy in units of electron mass
821  G4double p; // Electron momentum in units of electron mass
822  G4double f; // Spectral shape function value
823  for (G4int ptn = 0; ptn < npti; ptn++) {
824  // Calculate simple phase space spectrum
825  e = 1. + e0*(ptn+0.5)/100.;
826  p = std::sqrt(e*e - 1.);
827  f = p*e*(e0-e+1)*(e0-e+1);
828 
829  // Apply Fermi factor to get allowed shape
830  f *= corrections.FermiFunction(e);
831 
832  // Apply shape factor for forbidden transitions
833  f *= corrections.ShapeFactor(betaType, p, e0-e+1.);
834  pdf[ptn] = f;
835  }
836 
837  G4RandGeneral* aRandomEnergy = new G4RandGeneral( pdf, npti);
838  G4BetaMinusDecayChannel *aBetaMinusChannel = new
839  G4BetaMinusDecayChannel(GetVerboseLevel(), &theParentNucleus,
840  b, c*MeV, a*MeV, 0, FBeta, aRandomEnergy);
841  aBetaMinusChannel->SetICM(applyICM);
842  aBetaMinusChannel->SetARM(applyARM);
843  aBetaMinusChannel->SetHLThreshold(halflifethreshold);
844  theDecayTable->Insert(aBetaMinusChannel);
845  modeSumBR[1] += b;
846  delete[] pdf;
847  } // c > 0
848  } // if not first record
849  }
850  break;
851 
852  case BetaPlus:
853  {
854  if (modeFirstRecord[2]) {
855  modeFirstRecord[2] = false;
856  modeTotalBR[2] = b;
857  } else {
858  e0 = c*MeV/0.510999 - 2.;
859  // Need to test e0 for nuclei which have Q < 2Me in their
860  // data files (e.g. z67.a162)
861  if (e0 > 0.) {
862  G4BetaDecayCorrections corrections(1-Z, A);
863 
864  // array to store function shape
865  G4int npti = 100;
866  G4double* pdf = new G4double[npti];
867 
868  G4double e; // Total positron energy in units of electron mass
869  G4double p; // Positron momentum in units of electron mass
870  G4double f; // Spectral shape function value
871  for (G4int ptn = 0; ptn < npti; ptn++) {
872  // Calculate simple phase space spectrum
873  e = 1. + e0*(ptn+0.5)/100.;
874  p = std::sqrt(e*e - 1.);
875  f = p*e*(e0-e+1)*(e0-e+1);
876 
877  // Apply Fermi factor to get allowed shape
878  f *= corrections.FermiFunction(e);
879 
880  // Apply shape factor for forbidden transitions
881  f *= corrections.ShapeFactor(betaType, p, e0-e+1.);
882  pdf[ptn] = f;
883  }
884  G4RandGeneral* aRandomEnergy = new G4RandGeneral(pdf, npti);
885  G4BetaPlusDecayChannel* aBetaPlusChannel = new
887  &theParentNucleus, b,
888  (c-1.021998)*MeV, a*MeV, 0,
889  FBeta, aRandomEnergy);
890  aBetaPlusChannel->SetICM(applyICM);
891  aBetaPlusChannel->SetARM(applyARM);
892  aBetaPlusChannel->SetHLThreshold(halflifethreshold);
893  theDecayTable->Insert(aBetaPlusChannel);
894  modeSumBR[2] += b;
895  delete[] pdf;
896  } // if e0 > 0
897  } // if not first record
898  }
899  break;
900 
901  case KshellEC: // K-shell electron capture
902 
903  if (modeFirstRecord[3]) {
904  modeFirstRecord[3] = false;
905  modeTotalBR[3] = b;
906  } else {
907  G4KshellECDecayChannel* aKECChannel =
909  &theParentNucleus,
910  b, c*MeV, a*MeV);
911  aKECChannel->SetICM(applyICM);
912  aKECChannel->SetARM(applyARM);
913  aKECChannel->SetHLThreshold(halflifethreshold);
914  theDecayTable->Insert(aKECChannel);
915  modeSumBR[3] += b;
916  }
917  break;
918 
919  case LshellEC: // L-shell electron capture
920 
921  if (modeFirstRecord[4]) {
922  modeFirstRecord[4] = false;
923  modeTotalBR[4] = b;
924  } else {
925  G4LshellECDecayChannel *aLECChannel = new
926  G4LshellECDecayChannel (GetVerboseLevel(), &theParentNucleus,
927  b, c*MeV, a*MeV);
928  aLECChannel->SetICM(applyICM);
929  aLECChannel->SetARM(applyARM);
930  aLECChannel->SetHLThreshold(halflifethreshold);
931  theDecayTable->Insert(aLECChannel);
932  modeSumBR[4] += b;
933  }
934  break;
935 
936  case MshellEC: // M-shell electron capture
937  // In this implementation it is added to L-shell case
938  if (modeFirstRecord[5]) {
939  modeFirstRecord[5] = false;
940  modeTotalBR[5] = b;
941  } else {
942  G4MshellECDecayChannel* aMECChannel =
944  &theParentNucleus,
945  b, c*MeV, a*MeV);
946  aMECChannel->SetICM(applyICM);
947  aMECChannel->SetARM(applyARM);
948  aMECChannel->SetHLThreshold(halflifethreshold);
949  theDecayTable->Insert(aMECChannel);
950  modeSumBR[5] += b;
951  }
952  break;
953 
954  case Alpha:
955  if (modeFirstRecord[6]) {
956  modeFirstRecord[6] = false;
957  modeTotalBR[6] = b;
958  } else {
959  G4AlphaDecayChannel* anAlphaChannel =
961  &theParentNucleus,
962  b, c*MeV, a*MeV);
963  anAlphaChannel->SetICM(applyICM);
964  anAlphaChannel->SetARM(applyARM);
965  anAlphaChannel->SetHLThreshold(halflifethreshold);
966  theDecayTable->Insert(anAlphaChannel);
967  modeSumBR[6] += b;
968  }
969  break;
970  case SpFission:
971  //Still needed to be implemented
972  //G4cout<<"Sp fission channel"<<a<<'\t'<<b<<'\t'<<c<<std::endl;
973  break;
974  case RDM_ERROR:
975 
976  default:
977  G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_000",
978  FatalException, "Selected decay mode does not exist");
979  } // switch
980  } // if char == W
981  } // if char == P
982  } // if char != #
983  } // While
984 
985  // Go through the decay table and make sure that the branching ratios are
986  // correctly normalised.
987 
988  G4VDecayChannel* theChannel = 0;
989  G4NuclearDecayChannel* theNuclearDecayChannel = 0;
990  G4String mode = "";
991 
992  G4double theBR = 0.0;
993  for (G4int i = 0; i < theDecayTable->entries(); i++) {
994  theChannel = theDecayTable->GetDecayChannel(i);
995  theNuclearDecayChannel = static_cast<G4NuclearDecayChannel*>(theChannel);
996  theDecayMode = theNuclearDecayChannel->GetDecayMode();
997 
998  if (theDecayMode != IT) {
999  theBR = theChannel->GetBR();
1000  theChannel->SetBR(theBR*modeTotalBR[theDecayMode]/modeSumBR[theDecayMode]);
1001  }
1002  }
1003  } // if (DecaySchemeFile)
1004  DecaySchemeFile.close();
1005 
1006  if (!found && lvl > 0) {
1007  // Case where IT cascade for excited isotopes has no entries in RDM database
1008  // Decay mode is isomeric transition.
1009  G4ITDecayChannel* anITChannel = new G4ITDecayChannel
1010  (GetVerboseLevel(), (const G4Ions*) &theParentNucleus, 1);
1011  anITChannel->SetICM(applyICM);
1012  anITChannel->SetARM(applyARM);
1013  anITChannel->SetHLThreshold(halflifethreshold);
1014  theDecayTable->Insert(anITChannel);
1015  }
1016  if (!theDecayTable) {
1017 
1018  // There is no radioactive decay data for this nucleus. Return a null
1019  // decay table.
1020  G4cerr << "G4RadoactiveDecay::LoadDecayTable() : cannot find ion radioactive decay file "
1021  << G4endl;
1022  theDecayTable = 0;
1023  return theDecayTable;
1024  }
1025 
1026  if (theDecayTable && GetVerboseLevel() > 1) {
1027  G4cout << "G4RadioactiveDecay::LoadDecayTable()" << G4endl;
1028  theDecayTable->DumpInfo();
1029  }
1030 
1031 #ifdef G4MULTITHREADED
1032  (*master_dkmap)[key] = theDecayTable; // store in master library
1033 #endif
1034  return theDecayTable;
1035 }
1036 
1037 
1039 {
1040  if (Z < 1 || A < 2) G4cout << "Z and A not valid!" << G4endl;
1041 
1042  std::ifstream DecaySchemeFile(filename);
1043  if (DecaySchemeFile) {
1044  G4int ID_ion = A*1000 + Z;
1045  theUserRadioactiveDataFiles[ID_ion] = filename;
1046  theIsotopeTable->AddUserDecayDataFile(Z,A,filename);
1047  } else {
1048  G4cout << "The file " << filename << " does not exist!" << G4endl;
1049  }
1050 }
1051 
1052 
1053 void
1055  G4int theG, std::vector<G4double> theRates,
1056  std::vector<G4double> theTaos)
1057 {
1058  //fill the decay rate vector
1059  theDecayRate.SetZ(theZ);
1060  theDecayRate.SetA(theA);
1061  theDecayRate.SetE(theE);
1062  theDecayRate.SetGeneration(theG);
1063  theDecayRate.SetDecayRateC(theRates);
1064  theDecayRate.SetTaos(theTaos);
1065 }
1066 
1067 
1068 void
1070 {
1071  // 1) To calculate all the coefficiecies required to derive the
1072  // radioactivities for all progeny of theParentNucleus
1073  //
1074  // 2) Add the coefficiencies to the decay rate table vector
1075  //
1076 
1077  //
1078  // Create and initialise variables used in the method.
1079  //
1080  theDecayRateVector.clear();
1081 
1082  G4int nGeneration = 0;
1083  std::vector<G4double> rates;
1084  std::vector<G4double> taos;
1085 
1086  // start rate is -1.
1087  // Eq.4.26 of the Technical Note
1088  rates.push_back(-1.);
1089  //
1090  //
1091  G4int A = ((const G4Ions*)(&theParentNucleus))->GetAtomicMass();
1092  G4int Z = ((const G4Ions*)(&theParentNucleus))->GetAtomicNumber();
1093  G4double E = ((const G4Ions*)(&theParentNucleus))->GetExcitationEnergy();
1094  G4double tao = theParentNucleus.GetPDGLifeTime();
1095  if (tao < 0.) tao = 1e-100;
1096  taos.push_back(tao);
1097  G4int nEntry = 0;
1098 
1099  //fill the decay rate with the intial isotope data
1100  SetDecayRate(Z,A,E,nGeneration,rates,taos);
1101 
1102  // store the decay rate in decay rate vector
1103  theDecayRateVector.push_back(theDecayRate);
1104  nEntry++;
1105 
1106  // now start treating the sencondary generations..
1107 
1108  G4bool stable = false;
1109  G4int i;
1110  G4int j;
1111  G4VDecayChannel* theChannel = 0;
1112  G4NuclearDecayChannel* theNuclearDecayChannel = 0;
1113  G4ITDecayChannel* theITChannel = 0;
1114  G4BetaMinusDecayChannel *theBetaMinusChannel = 0;
1115  G4BetaPlusDecayChannel *theBetaPlusChannel = 0;
1116  G4AlphaDecayChannel *theAlphaChannel = 0;
1117  G4RadioactiveDecayMode theDecayMode;
1118  G4double theBR = 0.0;
1119  G4int AP = 0;
1120  G4int ZP = 0;
1121  G4int AD = 0;
1122  G4int ZD = 0;
1123  G4double EP = 0.;
1124  std::vector<G4double> TP;
1125  std::vector<G4double> RP;
1126  G4ParticleDefinition *theDaughterNucleus;
1127  G4double daughterExcitation;
1128  G4ParticleDefinition *aParentNucleus;
1129  G4IonTable* theIonTable;
1130  G4DecayTable *aTempDecayTable;
1131  G4double theRate;
1132  G4double TaoPlus;
1133  G4int nS = 0;
1134  G4int nT = nEntry;
1135  G4double brs[7];
1136  //
1137  theIonTable =
1139 
1140  while (!stable) {
1141  nGeneration++;
1142  for (j = nS; j < nT; j++) {
1143  ZP = theDecayRateVector[j].GetZ();
1144  AP = theDecayRateVector[j].GetA();
1145  EP = theDecayRateVector[j].GetE();
1146  RP = theDecayRateVector[j].GetDecayRateC();
1147  TP = theDecayRateVector[j].GetTaos();
1148  if (GetVerboseLevel() > 0) {
1149  G4cout << "G4RadioactiveDecay::AddDecayRateTable : daughters of ("
1150  << ZP << ", " << AP << ", " << EP
1151  << ") are being calculated, generation = " << nGeneration
1152  << G4endl;
1153  }
1154 
1155  aParentNucleus = theIonTable->GetIon(ZP,AP,EP);
1156  aTempDecayTable = GetDecayTable(aParentNucleus);
1157 
1158  G4DecayTable* theDecayTable = new G4DecayTable();
1159  for (G4int k = 0; k < 7; k++) brs[k] = 0.0;
1160 
1161  // Go through the decay table and to combine the same decay channels
1162  for (i = 0; i < aTempDecayTable->entries(); i++) {
1163  theChannel = aTempDecayTable->GetDecayChannel(i);
1164  theNuclearDecayChannel = static_cast<G4NuclearDecayChannel*>(theChannel);
1165  theDecayMode = theNuclearDecayChannel->GetDecayMode();
1166  daughterExcitation = theNuclearDecayChannel->GetDaughterExcitation();
1167  theDaughterNucleus = theNuclearDecayChannel->GetDaughterNucleus() ;
1168  AD = ((const G4Ions*)(theDaughterNucleus))->GetAtomicMass();
1169  ZD = ((const G4Ions*)(theDaughterNucleus))->GetAtomicNumber();
1170  G4NuclearLevelManager* levelManager =
1172  if (levelManager->NumberOfLevels() ) {
1173  const G4NuclearLevel* level =
1174  levelManager->NearestLevel (daughterExcitation);
1175 
1176  if (std::abs(daughterExcitation - level->Energy()) < levelTolerance) {
1177  // Level half-life is in ns and the threshold is set to 1 micros
1178  // by default, user can set it via the UI command
1179  if (level->HalfLife()*ns >= halflifethreshold){
1180  // save the metastable nucleus
1181  theDecayTable->Insert(theChannel);
1182  } else {
1183  brs[theDecayMode] += theChannel->GetBR();
1184  }
1185  } else {
1186  brs[theDecayMode] += theChannel->GetBR();
1187  }
1188  } else {
1189  brs[theDecayMode] += theChannel->GetBR();
1190  }
1191  }
1192  brs[2] = brs[2]+brs[3]+brs[4]+brs[5];
1193  brs[3] = brs[4] =brs[5] = 0.0;
1194  for (i= 0; i<7; i++){
1195  if (brs[i] > 0.) {
1196  switch ( i ) {
1197  case 0:
1198  // Decay mode is isomeric transition
1199  theITChannel = new G4ITDecayChannel(0,
1200  (const G4Ions*) aParentNucleus, brs[0]);
1201  theDecayTable->Insert(theITChannel);
1202  break;
1203 
1204  case 1:
1205  // Decay mode is beta-
1206  theBetaMinusChannel = new G4BetaMinusDecayChannel(0, aParentNucleus,
1207  brs[1], 0.*MeV, 0.*MeV, 1, false, 0);
1208  theDecayTable->Insert(theBetaMinusChannel);
1209  break;
1210 
1211  case 2:
1212  // Decay mode is beta+ + EC.
1213  theBetaPlusChannel = new G4BetaPlusDecayChannel(GetVerboseLevel(),
1214  aParentNucleus, brs[2], 0.*MeV, 0.*MeV, 1, false, 0);
1215  theDecayTable->Insert(theBetaPlusChannel);
1216  break;
1217 
1218  case 6:
1219  // Decay mode is alpha.
1220  theAlphaChannel = new G4AlphaDecayChannel(GetVerboseLevel(),
1221  aParentNucleus,
1222  brs[6], 0.*MeV, 0.*MeV);
1223  theDecayTable->Insert(theAlphaChannel);
1224  break;
1225 
1226  default:
1227  break;
1228  }
1229  }
1230  }
1231 
1232  // loop over all branches in theDecayTable
1233  //
1234  for (i = 0; i < theDecayTable->entries(); i++){
1235  theChannel = theDecayTable->GetDecayChannel(i);
1236  theNuclearDecayChannel = static_cast<G4NuclearDecayChannel*>(theChannel);
1237  theBR = theChannel->GetBR();
1238  theDaughterNucleus = theNuclearDecayChannel->GetDaughterNucleus();
1239  // First check if the decay of the original nucleus is an IT channel,
1240  // if true create a new groud-level nucleus
1241  if (theNuclearDecayChannel->GetDecayMode() == IT && nGeneration == 1) {
1242  A = ((const G4Ions*)(theDaughterNucleus))->GetAtomicMass();
1243  Z = ((const G4Ions*)(theDaughterNucleus))->GetAtomicNumber();
1244  theDaughterNucleus=theIonTable->GetIon(Z,A,0.);
1245  }
1246  if (IsApplicable(*theDaughterNucleus) && theBR &&
1247  aParentNucleus != theDaughterNucleus) {
1248  // need to make sure daughter has decay table
1249  aTempDecayTable = GetDecayTable(theDaughterNucleus);
1250 
1251  if (aTempDecayTable->entries() ) {
1252  A = ((const G4Ions*)(theDaughterNucleus))->GetAtomicMass();
1253  Z = ((const G4Ions*)(theDaughterNucleus))->GetAtomicNumber();
1254  E = ((const G4Ions*)(theDaughterNucleus))->GetExcitationEnergy();
1255 
1256  TaoPlus = theDaughterNucleus->GetPDGLifeTime();
1257  if (TaoPlus <= 0.) TaoPlus = 1e-100;
1258 
1259  // first set the taos, one simply need to add to the parent ones
1260  taos.clear();
1261  taos = TP;
1262  size_t k;
1263  //check that TaoPlus differs from other taos from at least 1.e5 relative difference
1264  //for (k = 0; k < TP.size(); k++){
1265  //if (std::abs((TaoPlus-TP[k])/TP[k])<1.e-5 ) TaoPlus=1.00001*TP[k];
1266  //}
1267  taos.push_back(TaoPlus);
1268  // now calculate the coefficiencies
1269  //
1270  // they are in two parts, first the less than n ones
1271  // Eq 4.24 of the TN
1272  rates.clear();
1273  long double ta1,ta2;
1274  ta2 = (long double)TaoPlus;
1275  for (k = 0; k < RP.size(); k++){
1276  ta1 = (long double)TP[k];
1277  if (ta1 == ta2) {
1278  theRate = 1.e100;
1279  } else {
1280  theRate = ta1/(ta1-ta2);
1281  }
1282  theRate = theRate * theBR * RP[k];
1283  rates.push_back(theRate);
1284  }
1285 
1286  // the sencond part: the n:n coefficiency
1287  // Eq 4.25 of the TN. Note Yn+1 is zero apart from Y1 which is -1
1288  // as treated at line 1013
1289  theRate = 0.;
1290  long double aRate, aRate1;
1291  aRate1 = 0.L;
1292  for (k = 0; k < RP.size(); k++){
1293  ta1 = (long double)TP[k];
1294  if (ta1 == ta2 ) {
1295  aRate = 1.e100;
1296  } else {
1297  aRate = ta2/(ta1-ta2);
1298  }
1299  aRate = aRate * (long double)(theBR * RP[k]);
1300  aRate1 += aRate;
1301  }
1302  theRate = -aRate1;
1303  rates.push_back(theRate);
1304  SetDecayRate (Z,A,E,nGeneration,rates,taos);
1305  theDecayRateVector.push_back(theDecayRate);
1306  nEntry++;
1307  }
1308  } // end of testing daughter nucleus
1309  } // end of i loop( the branches)
1310  // delete theDecayTable;
1311 
1312  } //end of for j loop
1313  nS = nT;
1314  nT = nEntry;
1315  if (nS == nT) stable = true;
1316  }
1317 
1318  // end of while loop
1319  // the calculation completed here
1320 
1321 
1322  // fill the first part of the decay rate table
1323  // which is the name of the original particle (isotope)
1324  theDecayRateTable.SetIonName(theParentNucleus.GetParticleName());
1325 
1326  // now fill the decay table with the newly completed decay rate vector
1327  theDecayRateTable.SetItsRates(theDecayRateVector);
1328 
1329  // finally add the decayratetable to the tablevector
1330  theDecayRateTableVector.push_back(theDecayRateTable);
1331 }
1332 
1334 // //
1335 // SetSourceTimeProfile //
1336 // read in the source time profile function (histogram) //
1337 // //
1339 
1341 {
1342  std::ifstream infile ( filename, std::ios::in );
1343  if (!infile) {
1345  ed << " Could not open file " << filename << G4endl;
1346  G4Exception("G4RadioactiveDecay::SetSourceTimeProfile()", "HAD_RDM_001",
1347  FatalException, ed);
1348  }
1349 
1350  G4double bin, flux;
1351  NSourceBin = -1;
1352  while (infile >> bin >> flux ) {
1353  NSourceBin++;
1354  if (NSourceBin > 99) {
1355  G4Exception("G4RadioactiveDecay::SetSourceTimeProfile()", "HAD_RDM_002",
1356  FatalException, "Input source time file too big (>100 rows)");
1357 
1358  } else {
1359  SBin[NSourceBin] = bin * s;
1360  SProfile[NSourceBin] = flux;
1361  }
1362  }
1364  infile.close();
1365 
1366 #ifdef G4VERBOSE
1367  if (GetVerboseLevel()>1)
1368  {G4cout <<" Source Timeprofile Nbin = " << NSourceBin <<G4endl;}
1369 #endif
1370 }
1371 
1373 // //
1374 // SetDecayBiasProfile //
1375 // read in the decay bias scheme function (histogram) //
1376 // //
1378 
1380 {
1381 
1382  std::ifstream infile(filename, std::ios::in);
1383  if (!infile) G4Exception("G4RadioactiveDecay::SetDecayBias()", "HAD_RDM_003",
1384  FatalException, "Unable to open bias data file" );
1385 
1386  G4double bin, flux;
1387  G4int dWindows = 0;
1388  G4int i ;
1389 
1390  theRadioactivityTables.clear();
1391 
1392  NDecayBin = -1;
1393  while (infile >> bin >> flux ) {
1394  NDecayBin++;
1395  if (NDecayBin > 99) {
1396  G4Exception("G4RadioactiveDecay::SetDecayBias()", "HAD_RDM_004",
1397  FatalException, "Input bias file too big (>100 rows)" );
1398  } else {
1399  DBin[NDecayBin] = bin * s;
1400  DProfile[NDecayBin] = flux;
1401  if (flux > 0.) {
1402  decayWindows[NDecayBin] = dWindows;
1403  dWindows++;
1404  G4RadioactivityTable *rTable = new G4RadioactivityTable() ;
1405  theRadioactivityTables.push_back(rTable);
1406  }
1407  }
1408  }
1409  for ( i = 1; i<= NDecayBin; i++) DProfile[i] += DProfile[i-1];
1410  for ( i = 0; i<= NDecayBin; i++) DProfile[i] /= DProfile[NDecayBin];
1411  // converted to accumulated probabilities
1412 
1414  infile.close();
1415 
1416 #ifdef G4VERBOSE
1417  if (GetVerboseLevel()>1)
1418  {G4cout <<" Decay Bias Profile Nbin = " << NDecayBin <<G4endl;}
1419 #endif
1420 }
1421 
1423 // //
1424 // DecayIt //
1425 // //
1427 
1430 {
1431  // Initialize G4ParticleChange object, get particle details and decay table
1432 
1433  fParticleChangeForRadDecay.Initialize(theTrack);
1434  const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
1435 
1436  G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
1437 
1438  // First check whether RDM applies to the current logical volume
1439  if (!isAllVolumesMode) {
1440  if (!std::binary_search(ValidVolumes.begin(), ValidVolumes.end(),
1441  theTrack.GetVolume()->GetLogicalVolume()->GetName())) {
1442 #ifdef G4VERBOSE
1443  if (GetVerboseLevel()>0) {
1444  G4cout <<"G4RadioactiveDecay::DecayIt : "
1445  << theTrack.GetVolume()->GetLogicalVolume()->GetName()
1446  << " is not selected for the RDM"<< G4endl;
1447  G4cout << " There are " << ValidVolumes.size() << " volumes" << G4endl;
1448  G4cout << " The Valid volumes are " << G4endl;
1449  for (size_t i = 0; i< ValidVolumes.size(); i++)
1450  G4cout << ValidVolumes[i] << G4endl;
1451  }
1452 #endif
1453  fParticleChangeForRadDecay.SetNumberOfSecondaries(0);
1454 
1455  // Kill the parent particle.
1456  fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill) ;
1457  fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(0.0);
1459  return &fParticleChangeForRadDecay;
1460  }
1461  }
1462 
1463  // Now check if particle is valid for RDM
1464  if (!(IsApplicable(*theParticleDef) ) ) {
1465  // Particle is not an ion or is outside the nucleuslimits for decay
1466 
1467 #ifdef G4VERBOSE
1468  if (GetVerboseLevel()>0) {
1469  G4cerr << "G4RadioactiveDecay::DecayIt : "
1470  << theParticleDef->GetParticleName()
1471  << " is not a valid nucleus for the RDM"<< G4endl;
1472  }
1473 #endif
1474  fParticleChangeForRadDecay.SetNumberOfSecondaries(0);
1475 
1476  // Kill the parent particle
1477  fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill) ;
1478  fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(0.0);
1480  return &fParticleChangeForRadDecay;
1481  }
1482 
1483  G4DecayTable* theDecayTable = GetDecayTable(theParticleDef);
1484 
1485  if (theDecayTable == 0 || theDecayTable->entries() == 0) {
1486  // No data in the decay table. Set particle change parameters
1487  // to indicate this.
1488 #ifdef G4VERBOSE
1489  if (GetVerboseLevel() > 0) {
1490  G4cerr <<"G4RadioactiveDecay::DecayIt : decay table not defined for ";
1491  G4cerr <<theParticleDef->GetParticleName() <<G4endl;
1492  }
1493 #endif
1494  fParticleChangeForRadDecay.SetNumberOfSecondaries(0);
1495 
1496  // Kill the parent particle.
1497  fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill) ;
1498  fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(0.0);
1500  return &fParticleChangeForRadDecay;
1501 
1502  } else {
1503  // Data found. Try to decay nucleus
1504  G4double energyDeposit = 0.0;
1505  G4double finalGlobalTime = theTrack.GetGlobalTime();
1506  G4double finalLocalTime = theTrack.GetLocalTime();
1507  G4int index;
1508  G4ThreeVector currentPosition;
1509  currentPosition = theTrack.GetPosition();
1510 
1511  // Check whether use Analogue or VR implementation
1512  if (AnalogueMC) {
1513 #ifdef G4VERBOSE
1514  if (GetVerboseLevel() > 0)
1515  G4cout <<"DecayIt: Analogue MC version " << G4endl;
1516 #endif
1517 
1518  G4DecayProducts* products = DoDecay(*theParticleDef);
1519 
1520  // Check if the product is the same as input and kill the track if
1521  // necessary to prevent infinite loop (11/05/10, F.Lei)
1522  if ( products->entries() == 1) {
1523  fParticleChangeForRadDecay.SetNumberOfSecondaries(0);
1524  fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill);
1525  fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(0.0);
1527  return &fParticleChangeForRadDecay;
1528  }
1529 
1530  // Get parent particle information and boost the decay products to the
1531  // laboratory frame based on this information.
1532 
1533  //The Parent Energy used for the boost should be the total energy of
1534  // the nucleus of the parent ion without the energy of the shell electrons
1535  // (correction for bug 1359 by L. Desorgher)
1536  G4double ParentEnergy = theParticle->GetKineticEnergy()
1537  + theParticle->GetParticleDefinition()->GetPDGMass();
1538  G4ThreeVector ParentDirection(theParticle->GetMomentumDirection());
1539 
1540  if (theTrack.GetTrackStatus() == fStopButAlive) {
1541  //this condition seems to be always True, further investigation is needed (L.Desorgher)
1542 
1543  // The particle is decayed at rest.
1544  // since the time is still for rest particle in G4 we need to add the
1545  // additional time lapsed between the particle come to rest and the
1546  // actual decay. This time is simply sampled with the mean-life of
1547  // the particle. But we need to protect the case PDGTime < 0.
1548  // (F.Lei 11/05/10)
1549  G4double temptime = -std::log( G4UniformRand())
1550  *theParticleDef->GetPDGLifeTime();
1551  if (temptime < 0.) temptime = 0.;
1552  finalGlobalTime += temptime;
1553  finalLocalTime += temptime;
1554  energyDeposit += theParticle->GetKineticEnergy();
1555  }
1556  products->Boost(ParentEnergy, ParentDirection);
1557 
1558  // Add products in theParticleChangeForRadDecay.
1559  G4int numberOfSecondaries = products->entries();
1560  fParticleChangeForRadDecay.SetNumberOfSecondaries(numberOfSecondaries);
1561 #ifdef G4VERBOSE
1562  if (GetVerboseLevel()>1) {
1563  G4cout <<"G4RadioactiveDecay::DecayIt : Decay vertex :";
1564  G4cout <<" Time: " <<finalGlobalTime/ns <<"[ns]";
1565  G4cout <<" X:" <<(theTrack.GetPosition()).x() /cm <<"[cm]";
1566  G4cout <<" Y:" <<(theTrack.GetPosition()).y() /cm <<"[cm]";
1567  G4cout <<" Z:" <<(theTrack.GetPosition()).z() /cm <<"[cm]";
1568  G4cout << G4endl;
1569  G4cout <<"G4Decay::DecayIt : decay products in Lab. Frame" <<G4endl;
1570  products->DumpInfo();
1571  products->IsChecked();
1572  }
1573 #endif
1574  for (index=0; index < numberOfSecondaries; index++) {
1575  G4Track* secondary = new G4Track(products->PopProducts(),
1576  finalGlobalTime, currentPosition);
1577  secondary->SetGoodForTrackingFlag();
1578  secondary->SetTouchableHandle(theTrack.GetTouchableHandle());
1579  fParticleChangeForRadDecay.AddSecondary(secondary);
1580  }
1581  delete products;
1582  // end of analogue MC algorithm
1583 
1584  } else {
1585  // Variance Reduction Method
1586 #ifdef G4VERBOSE
1587  if (GetVerboseLevel()>0)
1588  G4cout << "DecayIt: Variance Reduction version " << G4endl;
1589 #endif
1590  if (!IsRateTableReady(*theParticleDef)) {
1591  // if the decayrates are not ready, calculate them and
1592  // add to the rate table vector
1593  AddDecayRateTable(*theParticleDef);
1594  }
1595  //retrieve the rates
1596  GetDecayRateTable(*theParticleDef);
1597 
1598  // declare some of the variables required in the implementation
1599  G4ParticleDefinition* parentNucleus;
1600  G4IonTable* theIonTable;
1601  G4int PZ;
1602  G4int PA;
1603  G4double PE;
1604  G4String keyName;
1605  std::vector<G4double> PT;
1606  std::vector<G4double> PR;
1607  G4double taotime;
1608  long double decayRate;
1609 
1610  size_t i;
1611  size_t j;
1612  G4int numberOfSecondaries;
1613  G4int totalNumberOfSecondaries = 0;
1614  G4double currentTime = 0.;
1615  G4int ndecaych;
1616  G4DynamicParticle* asecondaryparticle;
1617  std::vector<G4DynamicParticle*> secondaryparticles;
1618  std::vector<G4double> pw;
1619  std::vector<G4double> ptime;
1620  pw.clear();
1621  ptime.clear();
1622 
1623  //now apply the nucleus splitting
1624  for (G4int n = 0; n < NSplit; n++) {
1625  // Get the decay time following the decay probability function
1626  // suppllied by user
1627  G4double theDecayTime = GetDecayTime();
1628  G4int nbin = GetDecayTimeBin(theDecayTime);
1629 
1630  // calculate the first part of the weight function
1631  G4double weight1 = 1.;
1632  if (nbin == 1) {
1633  weight1 = 1./DProfile[nbin-1]
1634  *(DBin[nbin]-DBin[nbin-1])/NSplit;
1635  } else if (nbin > 1) {
1636  weight1 = 1./(DProfile[nbin]-DProfile[nbin-2])
1637  *(DBin[nbin]-DBin[nbin-1])/NSplit;
1638  }
1639 
1640  // it should be calculated in seconds
1641  weight1 /= s ;
1642 
1643  // loop over all the possible secondaries of the nucleus
1644  // the first one is itself.
1645  for (i = 0; i < theDecayRateVector.size(); i++) {
1646  PZ = theDecayRateVector[i].GetZ();
1647  PA = theDecayRateVector[i].GetA();
1648  PE = theDecayRateVector[i].GetE();
1649  PT = theDecayRateVector[i].GetTaos();
1650  PR = theDecayRateVector[i].GetDecayRateC();
1651 
1652  // Calculate the decay rate of the isotope
1653  // decayRate is the radioactivity of isotope (PZ,PA,PE) at the
1654  // time 'theDecayTime'
1655  // it will be used to calculate the statistical weight of the
1656  // decay products of this isotope
1657 
1658  // G4cout <<"PA= "<< PA << " PZ= " << PZ << " PE= "<< PE <<G4endl;
1659  decayRate = 0.L;
1660  for (j = 0; j < PT.size(); j++) {
1661  taotime = GetTaoTime(theDecayTime,PT[j]);
1662  decayRate -= PR[j] * (long double)taotime;
1663  // Eq.4.23 of of the TN
1664  // note the negative here is required as the rate in the
1665  // equation is defined to be negative,
1666  // i.e. decay away, but we need positive value here.
1667 
1668  // G4cout << j << "\t"<< PT[j]/s <<"\t"<<PR[j]<< "\t"
1669  // << decayRate << G4endl;
1670  }
1671 
1672  // add the isotope to the radioactivity tables
1673  // G4cout <<theDecayTime/s <<"\t"<<nbin<<G4endl;
1674  // G4cout << theTrack.GetWeight() <<"\t"<<weight1<<"\t"<<decayRate<< G4endl;
1675  theRadioactivityTables[decayWindows[nbin-1]]->AddIsotope(PZ,PA,PE,weight1*decayRate,theTrack.GetWeight());
1676 
1677  // Now calculate the statistical weight
1678  // One needs to fold the source bias function with the decaytime
1679  // also need to include the track weight! (F.Lei, 28/10/10)
1680  G4double weight = weight1*decayRate*theTrack.GetWeight();
1681 
1682  // decay the isotope
1684  parentNucleus = theIonTable->GetIon(PZ,PA,PE);
1685 
1686  // Create a temprary products buffer.
1687  // Its contents to be transfered to the products at the end of the loop
1688  G4DecayProducts* tempprods = 0;
1689 
1690  // Decide whether to apply branching ratio bias or not
1691  if (BRBias) {
1692  G4DecayTable* decayTable = GetDecayTable(parentNucleus);
1693 
1694  ndecaych = G4int(decayTable->entries()*G4UniformRand());
1695  G4VDecayChannel* theDecayChannel = decayTable->GetDecayChannel(ndecaych);
1696  if (theDecayChannel == 0) {
1697  // Decay channel not found.
1698 #ifdef G4VERBOSE
1699  if (GetVerboseLevel()>0) {
1700  G4cerr << " G4RadioactiveDecay::DoIt : cannot determine decay channel ";
1701  G4cerr << " for this nucleus; decay as if no biasing active ";
1702  G4cerr << G4endl;
1703  decayTable ->DumpInfo();
1704  }
1705 #endif
1706  tempprods = DoDecay(*parentNucleus); // DHW 6 Dec 2010 - do decay as if no biasing
1707  // to avoid deref of temppprods = 0
1708  } else {
1709  // A decay channel has been identified, so execute the DecayIt.
1710  G4double tempmass = parentNucleus->GetPDGMass();
1711  tempprods = theDecayChannel->DecayIt(tempmass);
1712  weight *= (theDecayChannel->GetBR())*(decayTable->entries());
1713  }
1714  } else {
1715  tempprods = DoDecay(*parentNucleus);
1716  }
1717 
1718  // save the secondaries for buffers
1719  numberOfSecondaries = tempprods->entries();
1720  currentTime = finalGlobalTime + theDecayTime;
1721  for (index = 0; index < numberOfSecondaries; index++) {
1722  asecondaryparticle = tempprods->PopProducts();
1723  if (asecondaryparticle->GetDefinition()->GetBaryonNumber() < 5) {
1724  pw.push_back(weight);
1725  ptime.push_back(currentTime);
1726  secondaryparticles.push_back(asecondaryparticle);
1727  }
1728  }
1729  delete tempprods;
1730 
1731  } // end of i loop
1732  } // end of n loop
1733 
1734  // now deal with the secondaries in the two stl containers
1735  // and submmit them back to the tracking manager
1736  totalNumberOfSecondaries = pw.size();
1737  fParticleChangeForRadDecay.SetNumberOfSecondaries(totalNumberOfSecondaries);
1738  for (index=0; index < totalNumberOfSecondaries; index++) {
1739  G4Track* secondary = new G4Track(secondaryparticles[index],
1740  ptime[index], currentPosition);
1741  secondary->SetGoodForTrackingFlag();
1742  secondary->SetTouchableHandle(theTrack.GetTouchableHandle());
1743  secondary->SetWeight(pw[index]);
1744  fParticleChangeForRadDecay.AddSecondary(secondary);
1745  }
1746  // make sure the original track is set to stop and its kinematic energy collected
1747  //
1748  //theTrack.SetTrackStatus(fStopButAlive);
1749  //energyDeposit += theParticle->GetKineticEnergy();
1750 
1751  } // End of Variance Reduction
1752 
1753  // Kill the parent particle
1754  fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill) ;
1755  fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(energyDeposit);
1756  fParticleChangeForRadDecay.ProposeLocalTime(finalLocalTime);
1757  // Reset NumberOfInteractionLengthLeft.
1759 
1760  return &fParticleChangeForRadDecay ;
1761  }
1762 }
1763 
1764 
1767 {
1768  G4DecayProducts* products = 0;
1769  G4DecayTable* theDecayTable = GetDecayTable(&theParticleDef);
1770 
1771  // Choose a decay channel.
1772 #ifdef G4VERBOSE
1773  if (GetVerboseLevel() > 0) G4cout << "Select a channel..." << G4endl;
1774 #endif
1775 
1776  G4VDecayChannel* theDecayChannel = theDecayTable->SelectADecayChannel();
1777 
1778  if (theDecayChannel == 0) {
1779  // Decay channel not found.
1780  G4cerr << "G4RadioactiveDecay::DoIt : can not determine decay channel";
1781  G4cerr << G4endl;
1782  } else {
1783  // A decay channel has been identified, so execute the DecayIt.
1784 #ifdef G4VERBOSE
1785  if (GetVerboseLevel() > 1) {
1786  G4cerr << "G4RadioactiveDecay::DoIt : selected decay channel addr:";
1787  G4cerr << theDecayChannel << G4endl;
1788  }
1789 #endif
1790  G4double tempmass = theParticleDef.GetPDGMass();
1791  products = theDecayChannel->DecayIt(tempmass);
1792 
1793  // Apply directional bias if requested by user
1794  CollimateDecay(products);
1795  }
1796 
1797  return products;
1798 }
1799 
1800 
1801 // Apply directional bias for "visible" daughters (e+-, gamma, n, alpha)
1802 
1804  if (origin == forceDecayDirection) return; // No collimation requested
1805  if (180.*deg == forceDecayHalfAngle) return;
1806  if (0 == products || 0 == products->entries()) return;
1807 
1808 #ifdef G4VERBOSE
1809  if (GetVerboseLevel() > 0) G4cout << "Begin of CollimateDecay..." << G4endl;
1810 #endif
1811 
1812  // Particles suitable for directional biasing (for if-blocks below)
1816  static const G4ParticleDefinition* gamma = G4Gamma::Definition();
1817  static const G4ParticleDefinition* alpha = G4Alpha::Definition();
1818 
1819  G4ThreeVector newDirection; // Re-use to avoid memory churn
1820  for (G4int i=0; i<products->entries(); i++) {
1821  G4DynamicParticle* daughter = (*products)[i];
1822  const G4ParticleDefinition* daughterType =
1823  daughter->GetParticleDefinition();
1824  if (daughterType == electron || daughterType == positron ||
1825  daughterType == neutron || daughterType == gamma ||
1826  daughterType == alpha) CollimateDecayProduct(daughter);
1827  }
1828 }
1829 
1831 #ifdef G4VERBOSE
1832  if (GetVerboseLevel() > 1) {
1833  G4cout << "CollimateDecayProduct for daughter "
1834  << daughter->GetParticleDefinition()->GetParticleName() << G4endl;
1835  }
1836 #endif
1837 
1839  if (origin != collimate) daughter->SetMomentumDirection(collimate);
1840 }
1841 
1842 
1843 // Choose random direction within collimation cone
1844 
1846  if (origin == forceDecayDirection) return origin; // Don't do collimation
1847  if (forceDecayHalfAngle == 180.*deg) return origin;
1848 
1849  G4ThreeVector dir = forceDecayDirection;
1850 
1851  // Return direction offset by random throw
1852  if (forceDecayHalfAngle > 0.) {
1853  // Generate uniform direction around central axis
1854  G4double phi = 2.*pi*G4UniformRand();
1855  G4double cosMin = std::cos(forceDecayHalfAngle);
1856  G4double cosTheta = (1.-cosMin)*G4UniformRand() + cosMin; // [cosMin,1.)
1857 
1858  dir.setPhi(dir.phi()+phi);
1859  dir.setTheta(dir.theta()+std::acos(cosTheta));
1860  }
1861 
1862 #ifdef G4VERBOSE
1863  if (GetVerboseLevel()>1)
1864  G4cout << " ChooseCollimationDirection returns " << dir << G4endl;
1865 #endif
1866 
1867  return dir;
1868 }
void RegisterIsotopeTable(G4VIsotopeTable *table)
Definition: G4IonTable.cc:1440
G4double GetBR() const
void setPhi(double)
void SetBR(G4double value)
static G4LossTableManager * Instance()
G4double GetLocalTime() const
void SelectAVolume(const G4String aVolume)
std::map< G4String, G4DecayTable * > DecayTableMap
tuple a
Definition: test.py:11
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4String GetName() const
G4double GetKineticEnergy() const
Int_t index
Definition: macro.C:9
G4bool IsRateTableReady(const G4ParticleDefinition &)
const G4DynamicParticle * GetDynamicParticle() const
G4double HalfLife() const
const XML_Char * s
Definition: expat.h:262
void SetTaos(std::vector< G4double > value)
G4String strip(G4int strip_Type=trailing, char c=' ')
ifstream in
Definition: comparison.C:7
void SetDecayRateC(std::vector< G4double > value)
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4double Energy() const
static G4Alpha * Definition()
Definition: G4Alpha.cc:49
const G4ThreeVector & GetPosition() const
G4int GetZMax() const
void AddSecondary(G4Track *aSecondary)
G4TrackStatus GetTrackStatus() const
float bin[41]
Definition: plottest35.C:14
static G4Electron * Definition()
Definition: G4Electron.cc:49
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VDecayChannel * GetDecayChannel(G4int index) const
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4int GetVerboseLevel() const
void SetMomentumDirection(const G4ThreeVector &aDirection)
static G4NuclearLevelStore * GetInstance()
G4ParticleDefinition * GetDefinition() const
tuple x
Definition: test.py:50
void SetSourceTimeProfile(G4String filename)
subroutine sort(A, N)
Definition: dpm25nuc7.f:4670
G4bool IsApplicable(const G4ParticleDefinition &)
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
double weight
Definition: plottest35.C:25
TFile * file
static G4Positron * Definition()
Definition: G4Positron.cc:49
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
int G4int
Definition: G4Types.hh:78
G4int GetAMin() const
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:158
const G4String & GetParticleName() const
TFile f
Definition: plotHisto.C:6
void SetWeight(G4double aValue)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetItsRates(G4RadioactiveDecayRates arate)
G4double GetTotalMomentum() const
Double_t y
Definition: plot.C:279
G4int entries() const
tuple b
Definition: test.py:12
Char_t n[5]
virtual void Initialize(const G4Track &)
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
Definition: G4Ions.hh:51
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double >, std::vector< G4double >)
G4DecayTable * GetDecayTable(G4ParticleDefinition *)
Float_t Z
Definition: plot.C:39
G4double GetMass() const
G4int GetAMax() const
bool G4bool
Definition: G4Types.hh:79
void DumpInfo() const
#define G4RandGeneral
Definition: Randomize.hh:84
const G4ThreeVector & GetMomentumDirection() const
void DumpInfo() const
void CollimateDecayProduct(G4DynamicParticle *product)
G4RadioactiveDecayMode
static G4LogicalVolumeStore * GetInstance()
G4ThreeVector ChooseCollimationDirection() const
G4int GetDecayTimeBin(const G4double aDecayTime)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetParticleType() const
G4DecayProducts * DoDecay(G4ParticleDefinition &theParticleDef)
Definition: G4Step.hh:76
double phi() const
void setTheta(double)
const G4ParticleDefinition * GetParticleDefinition() const
int nanosecond
Definition: hepunit.py:82
G4double GetGlobalTime() const
G4double GetTaoTime(const G4double, const G4double)
const G4TouchableHandle & GetTouchableHandle() const
G4BetaDecayType
G4VDecayChannel * SelectADecayChannel()
Definition: G4DecayTable.cc:81
G4NuclearLevelManager * GetManager(G4int Z, G4int A)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double theta() const
G4ParticleDefinition * GetDaughterNucleus()
void Insert(G4VDecayChannel *aChannel)
Definition: G4DecayTable.cc:60
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4int G4Mutex
Definition: G4Threading.hh:156
G4LogicalVolume * GetLogicalVolume() const
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
void GetDecayRateTable(const G4ParticleDefinition &)
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * PopProducts()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4RadioactiveDecayMode GetDecayMode()
void AddDecayRateTable(const G4ParticleDefinition &)
#define DBL_MIN
Definition: templates.hh:75
static G4Neutron * Definition()
Definition: G4Neutron.cc:54
G4VPhysicalVolume * GetVolume() const
tuple z
Definition: test.py:28
G4double GetWeight() const
G4double GetPDGLifeTime() const
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay")
#define G4endl
Definition: G4ios.hh:61
G4bool IsChecked() const
G4VAtomDeexcitation * AtomDeexcitation()
G4DecayTable * LoadDecayTable(G4ParticleDefinition &theParentNucleus)
G4int entries() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double ShapeFactor(const G4BetaDecayType &, const G4double &p_e, const G4double &e_nu)
tuple c
Definition: test.py:13
G4ForceCondition
void SetHLThreshold(G4double hl)
void SetAnalogueMonteCarlo(G4bool r)
const G4NuclearLevel * NearestLevel(G4double energy, G4double eDiffMax=9999.*CLHEP::GeV) const
G4double FermiFunction(const G4double &W)
void DeselectAVolume(const G4String aVolume)
#define DBL_MAX
Definition: templates.hh:83
#define ns
Definition: xmlparse.cc:597
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
float c_light
Definition: hepunit.py:257
void CollimateDecay(G4DecayProducts *products)
void SetDecayBias(G4String filename)
G4int GetZMin() const
TDirectory * dir
Definition: macro.C:5
static G4Gamma * Definition()
Definition: G4Gamma.cc:49
void SetGoodForTrackingFlag(G4bool value=true)
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)
G4GLOB_DLL std::ostream G4cerr