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