Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Radioactivation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26////////////////////////////////////////////////////////////////////////////////
27// //
28// File: G4Radioactivation.cc //
29// Author: D.H. Wright (SLAC) //
30// Date: 29 August 2017 //
31// Description: activation process derived from the original //
32// G4RadioactiveDecay of F. Lei and P.R. Truscott in which //
33// biasing and activation calculations are separated from the //
34// unbiased decay chain calculation performed in the base //
35// class. //
36// //
37////////////////////////////////////////////////////////////////////////////////
38
39#include "G4Radioactivation.hh"
41
42#include "G4SystemOfUnits.hh"
43#include "G4DynamicParticle.hh"
44#include "G4DecayProducts.hh"
45#include "G4DecayTable.hh"
47#include "G4ITDecay.hh"
48#include "G4BetaDecayType.hh"
49#include "G4BetaMinusDecay.hh"
50#include "G4BetaPlusDecay.hh"
51#include "G4ECDecay.hh"
52#include "G4AlphaDecay.hh"
53#include "G4TritonDecay.hh"
54#include "G4ProtonDecay.hh"
55#include "G4NeutronDecay.hh"
56#include "G4SFDecay.hh"
57#include "G4VDecayChannel.hh"
58#include "G4NuclearDecay.hh"
60#include "G4Fragment.hh"
61#include "G4Ions.hh"
62#include "G4IonTable.hh"
63#include "G4BetaDecayType.hh"
64#include "Randomize.hh"
66#include "G4NuclearLevelData.hh"
68#include "G4LevelManager.hh"
69#include "G4ThreeVector.hh"
70#include "G4Electron.hh"
71#include "G4Positron.hh"
72#include "G4Neutron.hh"
73#include "G4Gamma.hh"
74#include "G4Alpha.hh"
75#include "G4Triton.hh"
76#include "G4Proton.hh"
77
81#include "G4LossTableManager.hh"
85
86#include <vector>
87#include <sstream>
88#include <algorithm>
89#include <fstream>
90
91using namespace CLHEP;
92
94 const G4double timeThresholdForRadioactiveDecays)
95 : G4RadioactiveDecay(processName, timeThresholdForRadioactiveDecays)
96{
97#ifdef G4VERBOSE
98 if (GetVerboseLevel() > 1) {
99 G4cout << "G4Radioactivation constructor: processName = " << processName
100 << G4endl;
101 }
102#endif
103
104 theRadioactivationMessenger = new G4RadioactivationMessenger(this);
105
106 // Apply default values.
107 NSourceBin = 1;
108 SBin[0] = 0.* s;
109 SBin[1] = 1.* s; // Convert to ns
110 SProfile[0] = 1.;
111 SProfile[1] = 0.;
112 NDecayBin = 1;
113 DBin[0] = 0. * s ;
114 DBin[1] = 1. * s;
115 DProfile[0] = 1.;
116 DProfile[1] = 0.;
117 decayWindows[0] = 0;
119 theRadioactivityTables.push_back(rTable);
120 NSplit = 1;
121 AnalogueMC = true;
122 BRBias = true;
123 halflifethreshold = 1000.*nanosecond;
124}
125
126void G4Radioactivation::ProcessDescription(std::ostream& outFile) const
127{
128 outFile << "The G4Radioactivation process performs radioactive decay of\n"
129 << "nuclides (G4GenericIon) in biased mode which includes nucleus\n"
130 << "duplication, branching ratio biasing, source time convolution\n"
131 << "and detector time convolution. It is designed for use in\n"
132 << "activation physics.\n"
133 << "The required half-lives and decay schemes are retrieved from\n"
134 << "the RadioactiveDecay database which was derived from ENSDF.\n";
135}
136
137
139{
140 delete theRadioactivationMessenger;
141}
142
143G4bool
145{
146 // Check whether the radioactive decay rates table for the ion has already
147 // been calculated.
148 G4String aParticleName = aParticle.GetParticleName();
149 for (std::size_t i = 0; i < theParentChainTable.size(); ++i) {
150 if (theParentChainTable[i].GetIonName() == aParticleName) return true;
151 }
152 return false;
153}
154
155void
157{
158 // Retrieve the decay rate table for the specified aParticle
159 G4String aParticleName = aParticle.GetParticleName();
160
161 for (std::size_t i = 0; i < theParentChainTable.size(); ++i) {
162 if (theParentChainTable[i].GetIonName() == aParticleName) {
163 theDecayRateVector = theParentChainTable[i].GetItsRates();
164 }
165 }
166#ifdef G4VERBOSE
167 if (GetVerboseLevel() > 1) {
168 G4cout << "The DecayRate Table for " << aParticleName << " is selected."
169 << G4endl;
170 }
171#endif
172}
173
174// ConvolveSourceTimeProfile performs the convolution of the source time profile
175// function with a single exponential characterized by a decay constant in the
176// decay chain. The time profile is treated as a step function so that the
177// convolution integral can be done bin-by-bin.
178// This implements Eq. 4.13 of DERA technical note, with SProfile[i] = F(t')
179
182{
183 G4double convolvedTime = 0.0;
184 G4int nbin;
185 if ( t > SBin[NSourceBin]) {
186 nbin = NSourceBin;
187 } else {
188 nbin = 0;
189
190 G4int loop = 0;
191 while (t > SBin[nbin]) { // Loop checking, 01.09.2015, D.Wright
192 loop++;
193 if (loop > 1000) {
194 G4Exception("G4Radioactivation::ConvolveSourceTimeProfile()",
195 "HAD_RDM_100", JustWarning, "While loop count exceeded");
196 break;
197 }
198 ++nbin;
199 }
200 --nbin;
201 }
202
203 // Use expm1 wherever possible to avoid large cancellation errors in
204 // 1 - exp(x) for small x
205 G4double earg = 0.0;
206 if (nbin > 0) {
207 for (G4int i = 0; i < nbin; ++i) {
208 earg = (SBin[i+1] - SBin[i])/tau;
209 if (earg < 100.) {
210 convolvedTime += SProfile[i] * std::exp((SBin[i] - t)/tau) *
211 std::expm1(earg);
212 } else {
213 convolvedTime += SProfile[i] *
214 (std::exp(-(t-SBin[i+1])/tau)-std::exp(-(t-SBin[i])/tau));
215 }
216 }
217 }
218 convolvedTime -= SProfile[nbin] * std::expm1((SBin[nbin] - t)/tau);
219 // tau divided out of final result to provide probability of decay in window
220
221 if (convolvedTime < 0.) {
222 G4cout << " Convolved time =: " << convolvedTime << " reset to zero! " << G4endl;
223 G4cout << " t = " << t << " tau = " << tau << G4endl;
224 G4cout << SBin[nbin] << " " << SBin[0] << G4endl;
225 convolvedTime = 0.;
226 }
227#ifdef G4VERBOSE
228 if (GetVerboseLevel() > 2)
229 G4cout << " Convolved time: " << convolvedTime << G4endl;
230#endif
231 return convolvedTime;
232}
233
234
235////////////////////////////////////////////////////////////////////////////////
236// //
237// GetDecayTime //
238// Randomly select a decay time for the decay process, following the //
239// supplied decay time bias scheme. //
240// //
241////////////////////////////////////////////////////////////////////////////////
242
244{
245 G4double decaytime = 0.;
246 G4double rand = G4UniformRand();
247 G4int i = 0;
248
249 G4int loop = 0;
250 while (DProfile[i] < rand) { /* Loop checking, 01.09.2015, D.Wright */
251 // Entries in DProfile[i] are all between 0 and 1 and arranged in inreaseing order
252 // Comparison with rand chooses which time bin to sample
253 ++i;
254 loop++;
255 if (loop > 100000) {
256 G4Exception("G4Radioactivation::GetDecayTime()", "HAD_RDM_100",
257 JustWarning, "While loop count exceeded");
258 break;
259 }
260 }
261
262 rand = G4UniformRand();
263 decaytime = DBin[i] + rand*(DBin[i+1]-DBin[i]);
264#ifdef G4VERBOSE
265 if (GetVerboseLevel() > 2)
266 G4cout <<" Decay time: " <<decaytime/s <<"[s]" <<G4endl;
267#endif
268 return decaytime;
269}
270
271
273{
274 G4int i = 0;
275
276 G4int loop = 0;
277 while (aDecayTime > DBin[i] ) { /* Loop checking, 01.09.2015, D.Wright */
278 ++i;
279 loop++;
280 if (loop > 100000) {
281 G4Exception("G4Radioactivation::GetDecayTimeBin()", "HAD_RDM_100",
282 JustWarning, "While loop count exceeded");
283 break;
284 }
285 }
286
287 return i;
288}
289
290////////////////////////////////////////////////////////////////////////////////
291// //
292// GetMeanLifeTime (required by the base class) //
293// //
294////////////////////////////////////////////////////////////////////////////////
295
298{
299 // For variance reduction time is set to 0 so as to force the particle
300 // to decay immediately.
301 // In analogue mode it returns the particle's mean-life.
302 G4double meanlife = 0.;
303 if (AnalogueMC) meanlife = G4RadioactiveDecay::GetMeanLifeTime(theTrack, fc);
304 return meanlife;
305}
306
307
308void
310 G4int theG, std::vector<G4double>& theCoefficients,
311 std::vector<G4double>& theTaos)
312// Why not make this a method of G4RadioactiveDecayRate? (e.g. SetParameters)
313{
314 //fill the decay rate vector
315 ratesToDaughter.SetZ(theZ);
316 ratesToDaughter.SetA(theA);
317 ratesToDaughter.SetE(theE);
318 ratesToDaughter.SetGeneration(theG);
319 ratesToDaughter.SetDecayRateC(theCoefficients);
320 ratesToDaughter.SetTaos(theTaos);
321}
322
323
325CalculateChainsFromParent(const G4ParticleDefinition& theParentNucleus)
326{
327 // Use extended Bateman equation to calculate the radioactivities of all
328 // progeny of theParentNucleus. The coefficients required to do this are
329 // calculated using the method of P. Truscott (Ph.D. thesis and
330 // DERA Technical Note DERA/CIS/CIS2/7/36/4/10) 11 January 2000.
331 // Coefficients are then added to the decay rate table vector
332
333 // Create and initialise variables used in the method.
334 theDecayRateVector.clear();
335
336 G4int nGeneration = 0;
337
338 std::vector<G4double> taos;
339
340 // Dimensionless A coefficients of Eqs. 4.24 and 4.25 of the TN
341 std::vector<G4double> Acoeffs;
342
343 // According to Eq. 4.26 the first coefficient (A_1:1) is -1
344 Acoeffs.push_back(-1.);
345
346 const G4Ions* ion = static_cast<const G4Ions*>(&theParentNucleus);
347 G4int A = ion->GetAtomicMass();
348 G4int Z = ion->GetAtomicNumber();
349 G4double E = ion->GetExcitationEnergy();
350 G4double tao = ion->GetPDGLifeTime();
351 if (tao < 0.) tao = 1e-100;
352 taos.push_back(tao);
353 G4int nEntry = 0;
354
355 // Fill the decay rate container (G4RadioactiveDecayRate) with the parent
356 // isotope data
357 SetDecayRate(Z,A,E,nGeneration,Acoeffs,taos); // Fill TP with parent lifetime
358
359 // store the decay rate in decay rate vector
360 theDecayRateVector.push_back(ratesToDaughter);
361 ++nEntry;
362
363 // Now start treating the secondary generations.
364 G4bool stable = false;
365 G4int j;
366 G4VDecayChannel* theChannel = 0;
367 G4NuclearDecay* theNuclearDecayChannel = 0;
368
369 G4ITDecay* theITChannel = 0;
370 G4BetaMinusDecay* theBetaMinusChannel = 0;
371 G4BetaPlusDecay* theBetaPlusChannel = 0;
372 G4AlphaDecay* theAlphaChannel = 0;
373 G4ProtonDecay* theProtonChannel = 0;
374 G4TritonDecay* theTritonChannel = 0;
375 G4NeutronDecay* theNeutronChannel = 0;
376 G4SFDecay* theFissionChannel = 0;
377
378 G4RadioactiveDecayMode theDecayMode;
379 G4double theBR = 0.0;
380 G4int AP = 0;
381 G4int ZP = 0;
382 G4int AD = 0;
383 G4int ZD = 0;
384 G4double EP = 0.;
385 std::vector<G4double> TP;
386 std::vector<G4double> RP; // A coefficients of the previous generation
387 G4ParticleDefinition *theDaughterNucleus;
388 G4double daughterExcitation;
389 G4double nearestEnergy = 0.0;
390 G4int nearestLevelIndex = 0;
391 G4ParticleDefinition *aParentNucleus;
392 G4IonTable* theIonTable;
393 G4DecayTable* parentDecayTable;
394 G4double theRate;
395 G4double TaoPlus;
396 G4int nS = 0; // Running index of first decay in a given generation
397 G4int nT = nEntry; // Total number of decays accumulated over entire history
398 const G4int nMode = G4RadioactiveDecayModeSize;
399 G4double brs[nMode];
400 //
402
403 G4int loop = 0;
404 while (!stable) { /* Loop checking, 01.09.2015, D.Wright */
405 loop++;
406 if (loop > 10000) {
407 G4Exception("G4Radioactivation::CalculateChainsFromParent()", "HAD_RDM_100",
408 JustWarning, "While loop count exceeded");
409 break;
410 }
411 nGeneration++;
412 for (j = nS; j < nT; ++j) {
413 // First time through, get data for parent nuclide
414 ZP = theDecayRateVector[j].GetZ();
415 AP = theDecayRateVector[j].GetA();
416 EP = theDecayRateVector[j].GetE();
417 RP = theDecayRateVector[j].GetDecayRateC();
418 TP = theDecayRateVector[j].GetTaos();
419 if (GetVerboseLevel() > 1) {
420 G4cout << "G4RadioactiveDecay::CalculateChainsFromParent: daughters of ("
421 << ZP << ", " << AP << ", " << EP
422 << ") are being calculated, generation = " << nGeneration
423 << G4endl;
424 }
425// G4cout << " Taus = " << G4endl;
426// for (G4int ii = 0; ii < TP.size(); ++ii) G4cout << TP[ii] << ", " ;
427// G4cout << G4endl;
428
429 aParentNucleus = theIonTable->GetIon(ZP,AP,EP);
430 parentDecayTable = GetDecayTable(aParentNucleus);
431 if (nullptr == parentDecayTable) { continue; }
432
433 G4DecayTable* summedDecayTable = new G4DecayTable();
434 // This instance of G4DecayTable is for accumulating BRs and decay
435 // channels. It will contain one decay channel per type of decay
436 // (alpha, beta, etc.); its branching ratio will be the sum of all
437 // branching ratios for that type of decay of the parent. If the
438 // halflife of a particular channel is longer than some threshold,
439 // that channel will be inserted specifically and its branching
440 // ratio will not be included in the above sums.
441 // This instance is not used to perform actual decays.
442
443 for (G4int k = 0; k < nMode; ++k) brs[k] = 0.0;
444
445 // Go through the decay table and sum all channels having the same decay mode
446 for (G4int i = 0; i < parentDecayTable->entries(); ++i) {
447 theChannel = parentDecayTable->GetDecayChannel(i);
448 theNuclearDecayChannel = static_cast<G4NuclearDecay*>(theChannel);
449 theDecayMode = theNuclearDecayChannel->GetDecayMode();
450 daughterExcitation = theNuclearDecayChannel->GetDaughterExcitation();
451 theDaughterNucleus = theNuclearDecayChannel->GetDaughterNucleus();
452 AD = ((const G4Ions*)(theDaughterNucleus))->GetAtomicMass();
453 ZD = ((const G4Ions*)(theDaughterNucleus))->GetAtomicNumber();
454 const G4LevelManager* levelManager =
456
457 // Check each nuclide to see if it is metastable (lifetime > 1 usec)
458 // If so, add it to the decay chain by inserting its decay channel in
459 // summedDecayTable. If not, just add its BR to sum for that decay mode.
460 if (levelManager->NumberOfTransitions() ) {
461 nearestEnergy = levelManager->NearestLevelEnergy(daughterExcitation);
462 if ((std::abs(daughterExcitation - nearestEnergy) < levelTolerance) && (std::abs(daughterExcitation - nearestEnergy) > DBL_EPSILON)) {
463 // Level half-life is in ns and the threshold is set to 1 micros
464 // by default, user can set it via the UI command
465 nearestLevelIndex = (G4int)levelManager->NearestLevelIndex(daughterExcitation);
466 if (levelManager->LifeTime(nearestLevelIndex)*ns >= halflifethreshold){
467 // save the metastable decay channel
468 summedDecayTable->Insert(theChannel);
469 } else {
470 brs[theDecayMode] += theChannel->GetBR();
471 }
472 } else {
473 brs[theDecayMode] += theChannel->GetBR();
474 }
475 } else {
476 brs[theDecayMode] += theChannel->GetBR();
477 }
478
479 } // Combine decay channels (loop i)
480
481 brs[BetaPlus] = brs[BetaPlus]+brs[KshellEC]+brs[LshellEC]+brs[MshellEC]+brs[NshellEC]; // Combine beta+ and EC
482 brs[KshellEC] = brs[LshellEC] = brs[MshellEC] = brs[NshellEC] = 0.0;
483 for (G4int i = 0; i < nMode; ++i) { // loop over decay modes
484 if (brs[i] > 0.) {
485 switch (i) {
486 case IT:
487 // Decay mode is isomeric transition
488 theITChannel = new G4ITDecay(aParentNucleus, brs[IT], 0.0, 0.0);
489
490 summedDecayTable->Insert(theITChannel);
491 break;
492
493 case BetaMinus:
494 // Decay mode is beta-
495 theBetaMinusChannel = new G4BetaMinusDecay(aParentNucleus, brs[BetaMinus],
496 0.*MeV, 0.*MeV,
498 summedDecayTable->Insert(theBetaMinusChannel);
499 break;
500
501 case BetaPlus:
502 // Decay mode is beta+ + EC.
503 theBetaPlusChannel = new G4BetaPlusDecay(aParentNucleus, brs[BetaPlus],
504 0.*MeV, 0.*MeV,
506 summedDecayTable->Insert(theBetaPlusChannel);
507 break;
508
509 case Alpha:
510 // Decay mode is alpha.
511 theAlphaChannel = new G4AlphaDecay(aParentNucleus, brs[Alpha], 0.*MeV,
512 0.*MeV, noFloat);
513 summedDecayTable->Insert(theAlphaChannel);
514 break;
515
516 case Proton:
517 // Decay mode is proton.
518 theProtonChannel = new G4ProtonDecay(aParentNucleus, brs[Proton], 0.*MeV,
519 0.*MeV, noFloat);
520 summedDecayTable->Insert(theProtonChannel);
521 break;
522
523 case Neutron:
524 // Decay mode is neutron.
525 theNeutronChannel = new G4NeutronDecay(aParentNucleus, brs[Neutron], 0.*MeV,
526 0.*MeV, noFloat);
527 summedDecayTable->Insert(theNeutronChannel);
528 break;
529
530 case SpFission:
531 // Decay mode is spontaneous fission
532 theFissionChannel = new G4SFDecay(aParentNucleus, brs[SpFission], 0.*MeV,
533 0.*MeV, noFloat);
534 summedDecayTable->Insert(theFissionChannel);
535 break;
536
537 case BDProton:
538 // Not yet implemented
539 break;
540
541 case BDNeutron:
542 // Not yet implemented
543 break;
544
545 case Beta2Minus:
546 // Not yet implemented
547 break;
548
549 case Beta2Plus:
550 // Not yet implemented
551 break;
552
553 case Proton2:
554 // Not yet implemented
555 break;
556
557 case Neutron2:
558 // Not yet implemented
559 break;
560
561 case Triton:
562 // Decay mode is Triton.
563 theTritonChannel = new G4TritonDecay(aParentNucleus, brs[Triton], 0.*MeV,
564 0.*MeV, noFloat);
565 summedDecayTable->Insert(theTritonChannel);
566 break;
567
568 default:
569 break;
570 }
571 }
572 }
573
574 // loop over all branches in summedDecayTable
575 //
576 for (G4int i = 0; i < summedDecayTable->entries(); ++i){
577 theChannel = summedDecayTable->GetDecayChannel(i);
578 theNuclearDecayChannel = static_cast<G4NuclearDecay*>(theChannel);
579 theBR = theChannel->GetBR();
580 theDaughterNucleus = theNuclearDecayChannel->GetDaughterNucleus();
581
582 // First check if the decay of the original nucleus is an IT channel,
583 // if true create a new ground-state nucleus
584 if (theNuclearDecayChannel->GetDecayMode() == IT && nGeneration == 1) {
585
586 A = ((const G4Ions*)(theDaughterNucleus))->GetAtomicMass();
587 Z = ((const G4Ions*)(theDaughterNucleus))->GetAtomicNumber();
588 theDaughterNucleus=theIonTable->GetIon(Z,A,0.);
589 }
590 if (IsApplicable(*theDaughterNucleus) && theBR > 0.0 &&
591 aParentNucleus != theDaughterNucleus) {
592 // need to make sure daughter has decay table
593 parentDecayTable = GetDecayTable(theDaughterNucleus);
594 if (nullptr != parentDecayTable && parentDecayTable->entries() > 0) {
595 A = ((const G4Ions*)(theDaughterNucleus))->GetAtomicMass();
596 Z = ((const G4Ions*)(theDaughterNucleus))->GetAtomicNumber();
597 E = ((const G4Ions*)(theDaughterNucleus))->GetExcitationEnergy();
598
599 TaoPlus = theDaughterNucleus->GetPDGLifeTime();
600 if (TaoPlus <= 0.) TaoPlus = 1e-100;
601
602 // first set the taos, one simply need to add to the parent ones
603 taos.clear();
604 taos = TP; // load lifetimes of all previous generations
605 std::size_t k;
606 //check that TaoPlus differs from other taos from at least 1.e5 relative difference
607 //for (k = 0; k < TP.size(); ++k){
608 //if (std::abs((TaoPlus-TP[k])/TP[k])<1.e-5 ) TaoPlus=1.00001*TP[k];
609 //}
610 taos.push_back(TaoPlus); // add daughter lifetime to list
611 // now calculate the coefficiencies
612 //
613 // they are in two parts, first the less than n ones
614 // Eq 4.24 of the TN
615 Acoeffs.clear();
616 long double ta1,ta2;
617 ta2 = (long double)TaoPlus;
618 for (k = 0; k < RP.size(); ++k){
619 ta1 = (long double)TP[k]; // loop over lifetimes of all previous generations
620 if (ta1 == ta2) {
621 theRate = 1.e100;
622 } else {
623 theRate = ta1/(ta1-ta2);
624 }
625 theRate = theRate * theBR * RP[k];
626 Acoeffs.push_back(theRate);
627 }
628
629 // the second part: the n:n coefficiency
630 // Eq 4.25 of the TN. Note Yn+1 is zero apart from Y1 which is -1
631 // as treated at line 1013
632 theRate = 0.;
633 long double aRate, aRate1;
634 aRate1 = 0.L;
635 for (k = 0; k < RP.size(); ++k){
636 ta1 = (long double)TP[k];
637 if (ta1 == ta2 ) {
638 aRate = 1.e100;
639 } else {
640 aRate = ta2/(ta1-ta2);
641 }
642 aRate = aRate * (long double)(theBR * RP[k]);
643 aRate1 += aRate;
644 }
645 theRate = -aRate1;
646 Acoeffs.push_back(theRate);
647 SetDecayRate (Z,A,E,nGeneration,Acoeffs,taos);
648 theDecayRateVector.push_back(ratesToDaughter);
649 nEntry++;
650 } // there are entries in the table
651 } // nuclide is OK to decay
652 } // end of loop (i) over decay table branches
653
654 delete summedDecayTable;
655
656 } // Getting contents of decay rate vector (end loop on j)
657 nS = nT;
658 nT = nEntry;
659 if (nS == nT) stable = true;
660 } // while nuclide is not stable
661
662 // end of while loop
663 // the calculation completed here
664
665
666 // fill the first part of the decay rate table
667 // which is the name of the original particle (isotope)
668 chainsFromParent.SetIonName(theParentNucleus.GetParticleName());
669
670 // now fill the decay table with the newly completed decay rate vector
671 chainsFromParent.SetItsRates(theDecayRateVector);
672
673 // finally add the decayratetable to the tablevector
674 theParentChainTable.push_back(chainsFromParent);
675}
676
677////////////////////////////////////////////////////////////////////////////////
678// //
679// SetSourceTimeProfile //
680// read in the source time profile function (histogram) //
681// //
682////////////////////////////////////////////////////////////////////////////////
683
685{
686 std::ifstream infile ( filename, std::ios::in );
687 if (!infile) {
689 ed << " Could not open file " << filename << G4endl;
690 G4Exception("G4Radioactivation::SetSourceTimeProfile()", "HAD_RDM_001",
691 FatalException, ed);
692 }
693
694 G4double bin, flux;
695 NSourceBin = -1;
696
697 G4int loop = 0;
698 while (infile >> bin >> flux) { /* Loop checking, 01.09.2015, D.Wright */
699 loop++;
700 if (loop > 10000) {
701 G4Exception("G4Radioactivation::SetSourceTimeProfile()", "HAD_RDM_100",
702 JustWarning, "While loop count exceeded");
703 break;
704 }
705
706 NSourceBin++;
707 if (NSourceBin > 99) {
708 G4Exception("G4Radioactivation::SetSourceTimeProfile()", "HAD_RDM_002",
709 FatalException, "Input source time file too big (>100 rows)");
710
711 } else {
712 SBin[NSourceBin] = bin * s; // Convert read-in time to ns
713 SProfile[NSourceBin] = flux; // Dimensionless
714 }
715 }
716
717 AnalogueMC = false;
718 infile.close();
719
720#ifdef G4VERBOSE
721 if (GetVerboseLevel() > 2)
722 G4cout <<" Source Timeprofile Nbin = " << NSourceBin <<G4endl;
723#endif
724}
725
726////////////////////////////////////////////////////////////////////////////////
727// //
728// SetDecayBiasProfile //
729// read in the decay bias scheme function (histogram) //
730// //
731////////////////////////////////////////////////////////////////////////////////
732
734{
735 std::ifstream infile(filename, std::ios::in);
736 if (!infile) G4Exception("G4Radioactivation::SetDecayBias()", "HAD_RDM_001",
737 FatalException, "Unable to open bias data file" );
738
739 G4double bin, flux;
740 G4int dWindows = 0;
741 G4int i ;
742
743 theRadioactivityTables.clear();
744
745 NDecayBin = -1;
746
747 G4int loop = 0;
748 while (infile >> bin >> flux ) { /* Loop checking, 01.09.2015, D.Wright */
749 NDecayBin++;
750 loop++;
751 if (loop > 10000) {
752 G4Exception("G4Radioactivation::SetDecayBias()", "HAD_RDM_100",
753 JustWarning, "While loop count exceeded");
754 break;
755 }
756
757 if (NDecayBin > 99) {
758 G4Exception("G4Radioactivation::SetDecayBias()", "HAD_RDM_002",
759 FatalException, "Input bias file too big (>100 rows)" );
760 } else {
761 DBin[NDecayBin] = bin * s; // Convert read-in time to ns
762 DProfile[NDecayBin] = flux; // Dimensionless
763 if (flux > 0.) {
764 decayWindows[NDecayBin] = dWindows;
765 dWindows++;
767 theRadioactivityTables.push_back(rTable);
768 }
769 }
770 }
771 for ( i = 1; i<= NDecayBin; ++i) DProfile[i] += DProfile[i-1]; // Cumulative flux vs i
772 for ( i = 0; i<= NDecayBin; ++i) DProfile[i] /= DProfile[NDecayBin];
773 // Normalize so entries increase from 0 to 1
774 // converted to accumulated probabilities
775
776 AnalogueMC = false;
777 infile.close();
778
779#ifdef G4VERBOSE
780 if (GetVerboseLevel() > 2)
781 G4cout <<" Decay Bias Profile Nbin = " << NDecayBin <<G4endl;
782#endif
783}
784
785////////////////////////////////////////////////////////////////////////////////
786// //
787// DecayIt //
788// //
789////////////////////////////////////////////////////////////////////////////////
790
793{
794 // Initialize G4ParticleChange object, get particle details and decay table
797 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
798 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
799
800 // First check whether RDM applies to the current logical volume
801 if (!isAllVolumesMode) {
802 if (!std::binary_search(ValidVolumes.begin(), ValidVolumes.end(),
803 theTrack.GetVolume()->GetLogicalVolume()->GetName())) {
804#ifdef G4VERBOSE
805 if (GetVerboseLevel()>0) {
806 G4cout <<"G4RadioactiveDecay::DecayIt : "
807 << theTrack.GetVolume()->GetLogicalVolume()->GetName()
808 << " is not selected for the RDM"<< G4endl;
809 G4cout << " There are " << ValidVolumes.size() << " volumes" << G4endl;
810 G4cout << " The Valid volumes are " << G4endl;
811 for (std::size_t i = 0; i< ValidVolumes.size(); ++i)
812 G4cout << ValidVolumes[i] << G4endl;
813 }
814#endif
816
817 // Kill the parent particle.
822 }
823 }
824
825 // Now check if particle is valid for RDM
826 if (!(IsApplicable(*theParticleDef) ) ) {
827 // Particle is not an ion or is outside the nucleuslimits for decay
828#ifdef G4VERBOSE
829 if (GetVerboseLevel() > 1) {
830 G4cout << "G4RadioactiveDecay::DecayIt : "
831 << theParticleDef->GetParticleName()
832 << " is not an ion or is outside (Z,A) limits set for the decay. "
833 << " Set particle change accordingly. "
834 << G4endl;
835 }
836#endif
838
839 // Kill the parent particle
844 }
845
846 G4DecayTable* theDecayTable = GetDecayTable(theParticleDef);
847
848 if (theDecayTable == nullptr || theDecayTable->entries() == 0) {
849 // No data in the decay table. Set particle change parameters
850 // to indicate this.
851#ifdef G4VERBOSE
852 if (GetVerboseLevel() > 1) {
853 G4cout << "G4RadioactiveDecay::DecayIt : "
854 << "decay table not defined for "
855 << theParticleDef->GetParticleName()
856 << ". Set particle change accordingly. "
857 << G4endl;
858 }
859#endif
861
862 // Kill the parent particle.
867
868 } else {
869 // Data found. Try to decay nucleus
870 if (AnalogueMC) {
871 G4RadioactiveDecay::DecayAnalog(theTrack, theDecayTable);
872
873 } else {
874 // Proceed with decay using variance reduction
875 G4double energyDeposit = 0.0;
876 G4double finalGlobalTime = theTrack.GetGlobalTime();
877 G4double finalLocalTime = theTrack.GetLocalTime();
878 G4int index;
879 G4ThreeVector currentPosition;
880 currentPosition = theTrack.GetPosition();
881
882 G4IonTable* theIonTable;
883 G4ParticleDefinition* parentNucleus;
884
885 // Get decay chains for the given nuclide
886 if (!IsRateTableReady(*theParticleDef))
887 CalculateChainsFromParent(*theParticleDef);
888 GetChainsFromParent(*theParticleDef);
889
890 // Declare some of the variables required in the implementation
891 G4int PZ;
892 G4int PA;
893 G4double PE;
894 G4String keyName;
895 std::vector<G4double> PT;
896 std::vector<G4double> PR;
897 G4double taotime;
898 long double decayRate;
899
900 std::size_t i;
901 G4int numberOfSecondaries;
902 G4int totalNumberOfSecondaries = 0;
903 G4double currentTime = 0.;
904 G4int ndecaych;
905 G4DynamicParticle* asecondaryparticle;
906 std::vector<G4DynamicParticle*> secondaryparticles;
907 std::vector<G4double> pw;
908 std::vector<G4double> ptime;
909 pw.clear();
910 ptime.clear();
911
912 // Now apply the nucleus splitting
913 for (G4int n = 0; n < NSplit; ++n) {
914 // Get the decay time following the decay probability function
915 // supplied by user
916 G4double theDecayTime = GetDecayTime();
917 G4int nbin = GetDecayTimeBin(theDecayTime);
918
919 // calculate the first part of the weight function
920 G4double weight1 = 1.;
921 if (nbin == 1) {
922 weight1 = 1./DProfile[nbin-1]
923 *(DBin[nbin]-DBin[nbin-1])/NSplit; // width of window in ns
924 } else if (nbin > 1) {
925 // Go from nbin to nbin-2 because flux entries in file alternate between 0 and 1
926 weight1 = 1./(DProfile[nbin]-DProfile[nbin-2])
927 *(DBin[nbin]-DBin[nbin-1])/NSplit;
928 // weight1 = (probability of choosing one of the bins)*(time width of bin)/NSplit
929 }
930 // it should be calculated in seconds
931 weight1 /= s ;
932
933 // loop over all the possible secondaries of the nucleus
934 // the first one is itself.
935 for (i = 0; i < theDecayRateVector.size(); ++i) {
936 PZ = theDecayRateVector[i].GetZ();
937 PA = theDecayRateVector[i].GetA();
938 PE = theDecayRateVector[i].GetE();
939 PT = theDecayRateVector[i].GetTaos();
940 PR = theDecayRateVector[i].GetDecayRateC();
941
942 // The array of arrays theDecayRateVector contains all possible decay
943 // chains of a given parent nucleus (ZP,AP,EP) to a given descendant
944 // nuclide (Z,A,E).
945 //
946 // theDecayRateVector[0] contains the decay parameters of the parent
947 // nucleus
948 // PZ = ZP
949 // PA = AP
950 // PE = EP
951 // PT[] = {TP}
952 // PR[] = {RP}
953 //
954 // theDecayRateVector[1] contains the decay of the parent to the first
955 // generation daughter (Z1,A1,E1).
956 // PZ = Z1
957 // PA = A1
958 // PE = E1
959 // PT[] = {TP, T1}
960 // PR[] = {RP, R1}
961 //
962 // theDecayRateVector[2] contains the decay of the parent to the first
963 // generation daughter (Z1,A1,E1) and the decay of the first
964 // generation daughter to the second generation daughter (Z2,A2,E2).
965 // PZ = Z2
966 // PA = A2
967 // PE = E2
968 // PT[] = {TP, T1, T2}
969 // PR[] = {RP, R1, R2}
970 //
971 // theDecayRateVector[3] may contain a branch chain
972 // PZ = Z2a
973 // PA = A2a
974 // PE = E2a
975 // PT[] = {TP, T1, T2a}
976 // PR[] = {RP, R1, R2a}
977 //
978 // and so on.
979
980 // Calculate the decay rate of the isotope. decayRate is the
981 // radioactivity of isotope (PZ,PA,PE) at 'theDecayTime'
982 // it will be used to calculate the statistical weight of the
983 // decay products of this isotope
984
985 // For each nuclide, calculate all the decay chains which can reach
986 // the parent nuclide
987 decayRate = 0.L;
988 for (G4int j = 0; j < G4int(PT.size() ); ++j) {
989 taotime = ConvolveSourceTimeProfile(theDecayTime,PT[j]);
990 decayRate -= PR[j] * (long double)taotime; // PRs are Acoeffs, taotime is inverse time
991 // Eq.4.23 of of the TN
992 // note the negative here is required as the rate in the
993 // equation is defined to be negative,
994 // i.e. decay away, but we need positive value here.
995
996 // G4cout << j << "\t"<< PT[j]/s << "\t" << PR[j] << "\t" << decayRate << G4endl;
997 }
998
999 // At this point any negative decay rates are probably small enough
1000 // (order 10**-30) that negative values are likely due to cancellation
1001 // errors. Set them to zero.
1002 if (decayRate < 0.0) decayRate = 0.0;
1003
1004 // G4cout <<theDecayTime/s <<"\t"<<nbin<<G4endl;
1005 // G4cout << theTrack.GetWeight() <<"\t"<<weight1<<"\t"<<decayRate<< G4endl;
1006
1007 // Add isotope to the radioactivity tables
1008 // One table for each observation time window specified in
1009 // SetDecayBias(G4String filename)
1010
1011 theRadioactivityTables[decayWindows[nbin-1]]
1012 ->AddIsotope(PZ,PA,PE,weight1*decayRate,theTrack.GetWeight());
1013
1014 // Now calculate the statistical weight
1015 // One needs to fold the source bias function with the decaytime
1016 // also need to include the track weight! (F.Lei, 28/10/10)
1017 G4double weight = weight1*decayRate*theTrack.GetWeight();
1018
1019 // decay the isotope
1020 theIonTable = (G4IonTable *)(G4ParticleTable::GetParticleTable()->GetIonTable());
1021 parentNucleus = theIonTable->GetIon(PZ,PA,PE);
1022
1023 // Create a temprary products buffer.
1024 // Its contents to be transfered to the products at the end of the loop
1025 G4DecayProducts* tempprods = nullptr;
1026
1027 // Decide whether to apply branching ratio bias or not
1028 if (BRBias) {
1029 G4DecayTable* decayTable = GetDecayTable(parentNucleus);
1030 G4VDecayChannel* theDecayChannel = nullptr;
1031 if (nullptr != decayTable) {
1032 ndecaych = G4int(decayTable->entries()*G4UniformRand());
1033 theDecayChannel = decayTable->GetDecayChannel(ndecaych);
1034 }
1035
1036 if (theDecayChannel == nullptr) {
1037 // Decay channel not found.
1038
1039 if (GetVerboseLevel() > 0) {
1040 G4cout << " G4RadioactiveDecay::DoIt : cannot determine decay channel ";
1041 G4cout << " for this nucleus; decay as if no biasing active. ";
1042 G4cout << G4endl;
1043 if (nullptr != decayTable) { decayTable ->DumpInfo(); }
1044 }
1045 // DHW 6 Dec 2010 - do decay as if no biasing to avoid deref of temppprods
1046 tempprods = DoDecay(*parentNucleus, theDecayTable);
1047 } else {
1048 // A decay channel has been identified, so execute the DecayIt.
1049 G4double tempmass = parentNucleus->GetPDGMass();
1050 tempprods = theDecayChannel->DecayIt(tempmass);
1051 weight *= (theDecayChannel->GetBR())*(decayTable->entries());
1052 }
1053 } else {
1054 tempprods = DoDecay(*parentNucleus, theDecayTable);
1055 }
1056
1057 // save the secondaries for buffers
1058 numberOfSecondaries = tempprods->entries();
1059 currentTime = finalGlobalTime + theDecayTime;
1060 for (index = 0; index < numberOfSecondaries; ++index) {
1061 asecondaryparticle = tempprods->PopProducts();
1062 if (asecondaryparticle->GetDefinition()->GetPDGStable() ) {
1063 pw.push_back(weight);
1064 ptime.push_back(currentTime);
1065 secondaryparticles.push_back(asecondaryparticle);
1066 }
1067 // Generate gammas and Xrays from excited nucleus, added by L.Desorgher
1068 else if (((const G4Ions*)(asecondaryparticle->GetDefinition()))
1069 ->GetExcitationEnergy() > 0. && weight > 0.) { //Compute the gamma
1070 G4ParticleDefinition* apartDef = asecondaryparticle->GetDefinition();
1071 AddDeexcitationSpectrumForBiasMode(apartDef,weight,currentTime,pw,
1072 ptime,secondaryparticles);
1073 }
1074 }
1075
1076 delete tempprods;
1077 } // end of i loop
1078 } // end of n loop
1079
1080 // now deal with the secondaries in the two stl containers
1081 // and submmit them back to the tracking manager
1082 totalNumberOfSecondaries = (G4int)pw.size();
1083 fParticleChangeForRadDecay.SetNumberOfSecondaries(totalNumberOfSecondaries);
1084 for (index=0; index < totalNumberOfSecondaries; ++index) {
1085 G4Track* secondary = new G4Track(secondaryparticles[index],
1086 ptime[index], currentPosition);
1087 secondary->SetGoodForTrackingFlag();
1088 secondary->SetTouchableHandle(theTrack.GetTouchableHandle());
1089 secondary->SetWeight(pw[index]);
1091 }
1092
1093 // Kill the parent particle
1097 // Reset NumberOfInteractionLengthLeft.
1099 } // end VR decay
1100
1102 } // end of data found branch
1103}
1104
1105
1106// Add gamma, X-ray, conversion and auger electrons for bias mode
1107void
1109 G4double weight,G4double currentTime,
1110 std::vector<double>& weights_v,
1111 std::vector<double>& times_v,
1112 std::vector<G4DynamicParticle*>& secondaries_v)
1113{
1114 G4double elevel=((const G4Ions*)(apartDef))->GetExcitationEnergy();
1115 G4double life_time=apartDef->GetPDGLifeTime();
1116 G4ITDecay* anITChannel = 0;
1117
1118 while (life_time < halflifethreshold && elevel > 0.) {
1119 decayIT->SetupDecay(apartDef);
1120 G4DecayProducts* pevap_products = decayIT->DecayIt(0.);
1121 G4int nb_pevapSecondaries = pevap_products->entries();
1122
1123 G4DynamicParticle* a_pevap_secondary = 0;
1124 G4ParticleDefinition* secDef = 0;
1125 for (G4int ind = 0; ind < nb_pevapSecondaries; ind++) {
1126 a_pevap_secondary= pevap_products->PopProducts();
1127 secDef = a_pevap_secondary->GetDefinition();
1128
1129 if (secDef->GetBaryonNumber() > 4) {
1130 elevel = ((const G4Ions*)(secDef))->GetExcitationEnergy();
1131 life_time = secDef->GetPDGLifeTime();
1132 apartDef = secDef;
1133 if (secDef->GetPDGStable() ) {
1134 weights_v.push_back(weight);
1135 times_v.push_back(currentTime);
1136 secondaries_v.push_back(a_pevap_secondary);
1137 }
1138 } else {
1139 weights_v.push_back(weight);
1140 times_v.push_back(currentTime);
1141 secondaries_v.push_back(a_pevap_secondary);
1142 }
1143 }
1144
1145 delete anITChannel;
1146 delete pevap_products;
1147 }
1148}
1149
@ allowed
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
#define noFloat
Definition G4Ions.hh:119
@ G4RadioactiveDecayModeSize
@ fStopAndKill
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4int entries() const
G4DynamicParticle * PopProducts()
G4VDecayChannel * GetDecayChannel(G4int index) const
void Insert(G4VDecayChannel *aChannel)
void DumpInfo() const
G4int entries() const
G4ParticleDefinition * GetDefinition() const
void SetupDecay(const G4ParticleDefinition *)
Definition G4ITDecay.cc:67
G4DecayProducts * DecayIt(G4double) override
Definition G4ITDecay.cc:74
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
G4double GetExcitationEnergy() const
Definition G4Ions.hh:149
std::size_t NearestLevelIndex(const G4double energy, const std::size_t index=0) const
std::size_t NumberOfTransitions() const
G4double LifeTime(const std::size_t i) const
G4double NearestLevelEnergy(const G4double energy, const std::size_t index=0) const
const G4String & GetName() const
G4RadioactiveDecayMode GetDecayMode() const
G4double GetDaughterExcitation() const
G4ParticleDefinition * GetDaughterNucleus()
const G4LevelManager * GetLevelManager(G4int Z, G4int A)
static G4NuclearLevelData * GetInstance()
void Initialize(const G4Track &) final
G4bool GetPDGStable() const
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4double GetPDGLifeTime() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
G4int GetDecayTimeBin(const G4double aDecayTime)
void AddDeexcitationSpectrumForBiasMode(G4ParticleDefinition *apartDef, G4double weight, G4double currenTime, std::vector< double > &weights_v, std::vector< double > &times_v, std::vector< G4DynamicParticle * > &secondaries_v)
G4bool IsRateTableReady(const G4ParticleDefinition &)
void SetDecayBias(const G4String &filename)
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition) override
G4double ConvolveSourceTimeProfile(const G4double, const G4double)
void CalculateChainsFromParent(const G4ParticleDefinition &)
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep) override
G4Radioactivation(const G4String &processName="Radioactivation", const G4double timeThresholdForRadioactiveDecays=-1.0)
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double > &, std::vector< G4double > &)
void ProcessDescription(std::ostream &outFile) const override
void SetSourceTimeProfile(const G4String &filename)
void GetChainsFromParent(const G4ParticleDefinition &)
void SetItsRates(G4RadioactiveDecayRates arate)
void SetTaos(std::vector< G4double > value)
void SetDecayRateC(std::vector< G4double > value)
G4bool IsApplicable(const G4ParticleDefinition &) override
std::vector< G4String > ValidVolumes
G4DecayProducts * DoDecay(const G4ParticleDefinition &, G4DecayTable *)
void DecayAnalog(const G4Track &theTrack, G4DecayTable *)
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
static const G4double levelTolerance
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition) override
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
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4LogicalVolume * GetLogicalVolume() const
G4int GetVerboseLevel() const
void ClearNumberOfInteractionLengthLeft()
#define DBL_EPSILON
Definition templates.hh:66