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