Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhotonEvaporation.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// GEANT4 class file
29//
30// CERN, Geneva, Switzerland
31//
32// File name: G4PhotonEvaporation
33//
34// Author: Vladimir Ivantchenko
35//
36// Creation date: 20 December 2011
37//
38//Modifications:
39//
40//
41// -------------------------------------------------------------------
42//
43
45
47#include "Randomize.hh"
48#include "G4Gamma.hh"
49#include "G4LorentzVector.hh"
50#include "G4FragmentVector.hh"
51#include "G4GammaTransition.hh"
52#include "G4Pow.hh"
55
56G4float G4PhotonEvaporation::GREnergy[] = {0.0f};
57G4float G4PhotonEvaporation::GRWidth[] = {0.0f};
58
59#ifdef G4MULTITHREADED
60G4Mutex G4PhotonEvaporation::PhotonEvaporationMutex = G4MUTEX_INITIALIZER;
61#endif
62
64 : fLevelManager(nullptr), fTransition(p), fPolarization(nullptr),
65 fVerbose(1), fPoints(0), vShellNumber(-1), fIndex(0),
66 fMaxLifeTime(DBL_MAX),
67 fICM(true), fRDM(false), fSampleTime(true),
68 fCorrelatedGamma(false), fIsomerFlag(false), isInitialised(false)
69{
70 //G4cout << "### New G4PhotonEvaporation() " << this << G4endl;
71 fNuclearLevelData = G4NuclearLevelData::GetInstance();
72 Tolerance = 20*CLHEP::eV;
73
74 if(!fTransition) { fTransition = new G4GammaTransition(); }
75
76 theA = theZ = fCode = 0;
77 fLevelEnergyMax = fStep = fExcEnergy = fProbability = 0.0;
78
79 for(G4int i=0; i<MAXDEPOINT; ++i) { fCummProbability[i] = 0.0; }
80 if(0.0f == GREnergy[1]) { InitialiseGRData(); }
81}
82
84{
85 delete fTransition;
86}
87
89{
90 if(isInitialised) { return; }
91 isInitialised = true;
92
93 G4DeexPrecoParameters* param = fNuclearLevelData->GetParameters();
94 Tolerance = param->GetMinExcitation();
95 fMaxLifeTime = param->GetMaxLifeTime();
96 fCorrelatedGamma = param->CorrelatedGamma();
97 fICM = param->GetInternalConversionFlag();
98 fIsomerFlag = param->IsomerProduction();
99 if(fRDM) { fIsomerFlag = true; }
100 fVerbose = param->GetVerbose();
101
102 fTransition->SetPolarizationFlag(fCorrelatedGamma);
103 fTransition->SetTwoJMAX(param->GetTwoJMAX());
104 fTransition->SetVerbose(fVerbose);
105 if(fVerbose > 1) {
106 G4cout << "### G4PhotonEvaporation is initialized " << this << G4endl;
107 }
108}
109
110void G4PhotonEvaporation::InitialiseGRData()
111{
112#ifdef G4MULTITHREADED
113 G4MUTEXLOCK(&G4PhotonEvaporation::PhotonEvaporationMutex);
114#endif
115 if(0.0f == GREnergy[1]) {
116 G4Pow* g4calc = G4Pow::GetInstance();
117 const G4float GRWfactor = 0.3f;
118 for (G4int A=1; A<MAXGRDATA; ++A) {
119 GREnergy[A] = (G4float)(40.3*CLHEP::MeV/g4calc->powZ(A,0.2));
120 GRWidth[A] = GRWfactor*GREnergy[A];
121 }
122 }
123#ifdef G4MULTITHREADED
124 G4MUTEXUNLOCK(&G4PhotonEvaporation::PhotonEvaporationMutex);
125#endif
126}
127
130{
131 if(!isInitialised) { Initialise(); }
132 fSampleTime = (fRDM) ? false : true;
133
134 // potentially external code may set initial polarization
135 // but only for radioactive decay nuclear polarization is considered
136 G4NuclearPolarizationStore* fNucPStore = nullptr;
137 if(fCorrelatedGamma && fRDM) {
139 if(nucleus->GetNuclearPolarization()) {
140 fNucPStore->RemoveMe(nucleus->GetNuclearPolarization());
141 delete nucleus->GetNuclearPolarization();
142 }
143 fPolarization = fNucPStore->FindOrBuild(nucleus->GetZ_asInt(),
144 nucleus->GetA_asInt(),
145 nucleus->GetExcitationEnergy());
146 nucleus->SetNuclearPolarization(fPolarization);
147 }
148 if(fVerbose > 2) {
149 G4cout << "G4PhotonEvaporation::EmittedFragment: "
150 << *nucleus << G4endl;
151 if(fPolarization) { G4cout << "NucPolar: " << fPolarization << G4endl; }
152 G4cout << " CorrGamma: " << fCorrelatedGamma << " RDM: " << fRDM
153 << " fPolarization: " << fPolarization << G4endl;
154 }
155 G4Fragment* gamma = GenerateGamma(nucleus);
156
157 // remove G4NuclearPolarizaton when reach ground state
158 if(fNucPStore && fPolarization && 0 == fIndex) {
159 if(fVerbose > 3) {
160 G4cout << "G4PhotonEvaporation::EmittedFragment: remove "
161 << fPolarization << G4endl;
162 }
163 fNucPStore->RemoveMe(fPolarization);
164 fPolarization = nullptr;
165 nucleus->SetNuclearPolarization(fPolarization);
166 }
167
168 if(fVerbose > 2) {
169 G4cout << "G4PhotonEvaporation::EmittedFragment: RDM= "
170 << fRDM << " done:" << G4endl;
171 if(gamma) { G4cout << *gamma << G4endl; }
172 G4cout << " Residual: " << *nucleus << G4endl;
173 }
174 return gamma;
175}
176
179{
180 if(fVerbose > 1) {
181 G4cout << "G4PhotonEvaporation::BreakItUp" << G4endl;
182 }
183 G4Fragment* aNucleus = new G4Fragment(nucleus);
184 G4FragmentVector* products = new G4FragmentVector();
185 BreakUpChain(products, aNucleus);
186 products->push_back(aNucleus);
187 return products;
188}
189
191 G4Fragment* nucleus)
192{
193 if(!isInitialised) { Initialise(); }
194 if(fVerbose > 1) {
195 G4cout << "G4PhotonEvaporation::BreakUpChain RDM= " << fRDM << " "
196 << *nucleus << G4endl;
197 }
198 G4Fragment* gamma = nullptr;
199 fSampleTime = (fRDM) ? false : true;
200
201 // start decay chain from unpolarized state
202 if(fCorrelatedGamma) {
203 fPolarization = new G4NuclearPolarization(nucleus->GetZ_asInt(),
204 nucleus->GetA_asInt(),
205 nucleus->GetExcitationEnergy());
206 nucleus->SetNuclearPolarization(fPolarization);
207 }
208
209 do {
210 gamma = GenerateGamma(nucleus);
211 if(gamma) {
212 products->push_back(gamma);
213 if(fVerbose > 2) {
214 G4cout << "G4PhotonEvaporation::BreakUpChain: "
215 << *gamma << G4endl;
216 G4cout << " Residual: " << *nucleus << G4endl;
217 }
218 // for next decays in the chain always sample time
219 fSampleTime = true;
220 }
221 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
222 } while(gamma);
223
224 // clear nuclear polarization end of chain
225 if(fPolarization) {
226 delete fPolarization;
227 fPolarization = nullptr;
228 nucleus->SetNuclearPolarization(fPolarization);
229 }
230 return false;
231}
232
235{
236 if(!isInitialised) { Initialise(); }
237 fProbability = 0.0;
238 fExcEnergy = nucleus->GetExcitationEnergy();
239 G4int Z = nucleus->GetZ_asInt();
240 G4int A = nucleus->GetA_asInt();
241 fCode = 1000*Z + A;
242 if(fVerbose > 2) {
243 G4cout << "G4PhotonEvaporation::GetEmissionProbability: Z="
244 << Z << " A=" << A << " Eexc(MeV)= " << fExcEnergy << G4endl;
245 }
246 // ignore gamma de-excitation for exotic fragments
247 // and for very low excitations
248 if(0 >= Z || 1 >= A || Z == A || Tolerance >= fExcEnergy)
249 { return fProbability; }
250
251 // ignore gamma de-excitation for highly excited levels
252 if(A >= MAXGRDATA) { A = MAXGRDATA-1; }
253 //G4cout<<" GREnergy= "<< GREnergy[A]<<" GRWidth= "<<GRWidth[A]<<G4endl;
254
255 static const G4float GREfactor = 5.0f;
256 if(fExcEnergy >= (G4double)(GREfactor*GRWidth[A] + GREnergy[A])) {
257 return fProbability;
258 }
259 // probability computed assuming continium transitions
260 // VI: continium transition are limited only to final states
261 // below Fermi energy (this approach needs further evaluation)
262 G4double emax = std::max(0.0, nucleus->ComputeGroundStateMass(Z, A-1)
263 + CLHEP::neutron_mass_c2 - nucleus->GetGroundStateMass());
264
265 // max energy level for continues transition
266 emax = std::min(emax, fExcEnergy);
267 const G4double eexcfac = 0.99;
268 if(0.0 == emax || fExcEnergy*eexcfac <= emax) { emax = fExcEnergy*eexcfac; }
269
270 fStep = emax;
271 const G4double MaxDeltaEnergy = CLHEP::MeV;
272 fPoints = std::min((G4int)(fStep/MaxDeltaEnergy) + 2, MAXDEPOINT);
273 fStep /= ((G4double)(fPoints - 1));
274 if(fVerbose > 2) {
275 G4cout << "Emax= " << emax << " Npoints= " << fPoints
276 << " Eex= " << fExcEnergy << G4endl;
277 }
278 // integrate probabilities
279 G4double eres = (G4double)GREnergy[A];
280 G4double wres = (G4double)GRWidth[A];
281 G4double eres2= eres*eres;
282 G4double wres2= wres*wres;
283 G4double levelDensity = fNuclearLevelData->GetLevelDensity(Z,A,fExcEnergy);
284 G4double xsqr = std::sqrt(levelDensity*fExcEnergy);
285
286 G4double egam = fExcEnergy;
287 G4double gammaE2 = egam*egam;
288 G4double gammaR2 = gammaE2*wres2;
289 G4double egdp2 = gammaE2 - eres2;
290
291 G4double p0 = G4Exp(-2.0*xsqr)*gammaR2*gammaE2/(egdp2*egdp2 + gammaR2);
292 G4double p1(0.0);
293
294 for(G4int i=1; i<fPoints; ++i) {
295 egam -= fStep;
296 gammaE2 = egam*egam;
297 gammaR2 = gammaE2*wres2;
298 egdp2 = gammaE2 - eres2;
299 p1 = G4Exp(2.0*(std::sqrt(levelDensity*std::abs(fExcEnergy - egam)) - xsqr))
300 *gammaR2*gammaE2/(egdp2*egdp2 + gammaR2);
301 fProbability += (p1 + p0);
302 fCummProbability[i] = fProbability;
303 if(fVerbose > 3) {
304 G4cout << "Egamma= " << egam << " Eex= " << fExcEnergy
305 << " p0= " << p0 << " p1= " << p1 << " sum= "
306 << fCummProbability[i] <<G4endl;
307 }
308 p0 = p1;
309 }
310
311 static const G4double NormC = 1.25*CLHEP::millibarn
312 /(CLHEP::pi2*CLHEP::hbarc*CLHEP::hbarc);
313 fProbability *= fStep*NormC*A;
314 if(fVerbose > 1) { G4cout << "prob= " << fProbability << G4endl; }
315 return fProbability;
316}
317
320{
321 G4double E = energy;
322 InitialiseLevelManager(Z, A);
323 if(fLevelManager) {
324 E = fLevelManager->NearestLevelEnergy(energy, fIndex);
325 if(E > fLevelEnergyMax + Tolerance) { E = energy; }
326 }
327 return E;
328}
329
331{
332 InitialiseLevelManager(Z, A);
333 return fLevelEnergyMax;
334}
335
337G4PhotonEvaporation::GenerateGamma(G4Fragment* nucleus)
338{
339 if(!isInitialised) { Initialise(); }
340 G4Fragment* result = nullptr;
341 G4double eexc = nucleus->GetExcitationEnergy();
342 if(eexc <= Tolerance) { return result; }
343
344 InitialiseLevelManager(nucleus->GetZ_asInt(), nucleus->GetA_asInt());
345
346 G4double time = nucleus->GetCreationTime();
347
348 G4double efinal = 0.0;
349 G4double ratio = 0.0;
350 vShellNumber = -1;
351 G4int JP1 = 0;
352 G4int JP2 = 0;
353 G4int multiP = 0;
354 G4bool isGamma = true;
355 G4bool isDiscrete = false;
356
357 const G4NucLevel* level = nullptr;
358 size_t ntrans = 0;
359
360 if(fVerbose > 2) {
361 G4cout << "GenerateGamma: " << " Eex= " << eexc
362 << " Eexmax= " << fLevelEnergyMax << G4endl;
363 }
364 // initial discrete state
365 if(fLevelManager && eexc <= fLevelEnergyMax + Tolerance) {
366 fIndex = fLevelManager->NearestLevelIndex(eexc, fIndex);
367 // initial state below 1st level
368 if(0 == fIndex && eexc >= Tolerance
369 && fLevelManager->NumberOfTransitions() > 0) { fIndex = 1; }
370 isDiscrete = true;
371 if(fVerbose > 2) {
372 G4cout << " index= " << fIndex
373 << " lTime= " << fLevelManager->LifeTime(fIndex) << G4endl;
374 }
375 if(0 < fIndex) {
376 // for discrete transition
377 level = fLevelManager->GetLevel(fIndex);
378 if(level) {
379 ntrans = level->NumberOfTransitions();
380 JP1 = fLevelManager->SpinTwo(fIndex);
381 if(fVerbose > 2) {
382 G4cout << " ntrans= " << ntrans << " JP= " << JP1
383 << " RDM: " << fRDM << G4endl;
384 }
385 if(0 == ntrans && fLevelManager->FloatingLevel(fIndex) > 0) {
386 --fIndex;
387 level = fLevelManager->GetLevel(fIndex);
388 ntrans = level->NumberOfTransitions();
389 JP1 = fLevelManager->SpinTwo(fIndex);
390 }
391 }
392 }
393 }
394 if(fVerbose > 2) {
395 G4int prec = G4cout.precision(4);
396 G4cout << "GenerateGamma: Z= " << nucleus->GetZ_asInt()
397 << " A= " << nucleus->GetA_asInt()
398 << " Exc= " << eexc << " Emax= "
399 << fLevelEnergyMax << " idx= " << fIndex
400 << " fCode= " << fCode << " fPoints= " << fPoints
401 << " Ntr= " << ntrans << " discrete: " << isDiscrete
402 << " fProb= " << fProbability << G4endl;
403 G4cout.precision(prec);
404 }
405
406 // continues part
407 if(!isDiscrete) {
408 // we compare current excitation versus value used for probability
409 // computation and also Z and A used for probability computation
410 if(fCode != 1000*theZ + theA || eexc != fExcEnergy) {
411 GetEmissionProbability(nucleus);
412 }
413 if(fProbability == 0.0) {
414 fPoints = 1;
415 efinal = 0.0;
416 } else {
417 G4double y = fCummProbability[fPoints-1]*G4UniformRand();
418 for(G4int i=1; i<fPoints; ++i) {
419 if(fVerbose > 3) {
420 G4cout << "y= " << y << " cummProb= " << fCummProbability[i]
421 << " fPoints= " << fPoints << " fStep= " << fStep << G4endl;
422 }
423 if(y <= fCummProbability[i]) {
424 efinal = fStep*((i - 1) + (y - fCummProbability[i-1])
425 /(fCummProbability[i] - fCummProbability[i-1]));
426 break;
427 }
428 }
429 }
430 // final discrete level
431 if(fVerbose > 2) {
432 G4cout << "Continues proposes Efinal= " << efinal << G4endl;
433 }
434 if(fLevelManager) {
435 if(efinal < fLevelEnergyMax) {
436 fIndex = fLevelManager->NearestLevelIndex(efinal, fIndex);
437 efinal = fLevelManager->LevelEnergy(fIndex);
438 // protection - take level below
439 if(efinal >= eexc && 0 < fIndex) {
440 --fIndex;
441 efinal = fLevelManager->LevelEnergy(fIndex);
442 }
443 nucleus->SetFloatingLevelNumber(fLevelManager->FloatingLevel(fIndex));
444
445 // not allowed to have final energy above max energy
446 // if G4LevelManager exist
447 } else {
448 efinal = fLevelEnergyMax;
449 fIndex = fLevelManager->NearestLevelIndex(efinal, fIndex);
450 }
451 }
452 if(fVerbose > 2) {
453 G4cout << "Continues emission efinal(MeV)= " << efinal << G4endl;
454 }
455 //discrete part ground state
456 } else if(0 == fIndex) {
457 return result;
458
459 //discrete part
460 } else {
461
462 if(fVerbose > 2) {
463 G4cout << "Discrete emission from level Index= " << fIndex
464 << " Elevel= " << fLevelManager->LevelEnergy(fIndex)
465 << " Ltime= " << fLevelManager->LifeTime(fIndex)
466 << " LtimeMax= " << fMaxLifeTime
467 << " RDM= " << fRDM << " ICM= " << fICM << G4endl;
468 }
469
470 // stable fragment has life time -1 or above the limit
471 // if is called from the radioactive decay the life time is not checked
472 G4double ltime = fLevelManager->LifeTime(fIndex);
473 if(ltime < 0.0 || (!fRDM && ltime > fMaxLifeTime)) { return result; }
474
475 // no transitions: force transition to the closest level
476 if(0 == ntrans) {
477 G4int ii = fIndex - 1;
478 for(; ii>0; --ii) {
479 const G4NucLevel* fl = fLevelManager->GetLevel(ii);
480 if(fl && 0 < fl->NumberOfTransitions()) { break; }
481 }
482 fIndex = ii;
483 } else {
484 size_t idx = 0;
485 if(1 < ntrans) {
486 idx = level->SampleGammaTransition(G4UniformRand());
487 }
488 if(fVerbose > 2) {
489 G4cout << "Ntrans= " << ntrans << " idx= " << idx
490 << " ICM= " << fICM << " JP1= " << JP1 << G4endl;
491 }
492 G4double prob = (G4double)level->GammaProbability(idx);
493 // prob = 0 means that there is only internal conversion
494 if(fICM && prob < 1.0) {
495 G4double rndm = G4UniformRand();
496 if(rndm > prob) {
497 isGamma = false;
498 rndm = (rndm - prob)/(1.0 - prob);
499 vShellNumber = level->SampleShell(idx, rndm);
500 }
501 }
502 // it is discrete transition with possible gamma correlation
503 ratio = level->MultipolarityRatio(idx);
504 multiP = level->TransitionType(idx);
505 fIndex = level->FinalExcitationIndex(idx);
506 }
507 JP2 = fLevelManager->SpinTwo(fIndex);
508
509 // final energy and time
510 efinal = fLevelManager->LevelEnergy(fIndex);
511 // time is sampled if decay not prompt and this class called not
512 // from radioactive decay and isomer production is enabled
513 if(fSampleTime && fIsomerFlag && ltime > 0.0) {
514 time -= ltime*G4Log(G4UniformRand());
515 }
516 nucleus->SetFloatingLevelNumber(fLevelManager->FloatingLevel(fIndex));
517 }
518 // protection for floating levels
519 if(std::abs(efinal - eexc) <= Tolerance) { return result; }
520
521 result = fTransition->SampleTransition(nucleus, efinal, ratio, JP1,
522 JP2, multiP, vShellNumber,
523 isDiscrete, isGamma);
524 if(result) { result->SetCreationTime(time); }
525
526 // updated residual nucleus
527 nucleus->SetCreationTime(time);
528 nucleus->SetSpin(0.5*JP2);
529 if(fPolarization) { fPolarization->SetExcitationEnergy(efinal); }
530
531 // ignore the floating levels with zero energy and create ground state
532 if(efinal == 0.0 && fIndex > 0) {
533 fIndex = 0;
534 nucleus->SetFloatingLevelNumber(fLevelManager->FloatingLevel(fIndex));
535 }
536
537 if(fVerbose > 2) {
538 G4cout << "Final level E= " << efinal << " time= " << time
539 << " idxFinal= " << fIndex << " isDiscrete: " << isDiscrete
540 << " isGamma: " << isGamma << " multiP= " << multiP
541 << " shell= " << vShellNumber
542 << " JP1= " << JP1 << " JP2= " << JP2 << G4endl;
543 }
544 return result;
545}
546
548{
549 if(p != fTransition) {
550 delete fTransition;
551 fTransition = p;
552 }
553}
554
556{
557 fICM = val;
558}
559
561{
562 fRDM = val;
563}
564
double A(double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
G4double G4Log(G4double x)
Definition: G4Log.hh:226
const G4int MAXDEPOINT
const G4int MAXGRDATA
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
float G4float
Definition: G4Types.hh:84
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetMinExcitation() const
G4bool GetInternalConversionFlag() const
void SetFloatingLevelNumber(G4int value)
Definition: G4Fragment.hh:424
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:280
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:455
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:275
G4double GetCreationTime() const
Definition: G4Fragment.hh:440
void SetNuclearPolarization(G4NuclearPolarization *)
Definition: G4Fragment.hh:460
G4int GetZ_asInt() const
Definition: G4Fragment.hh:263
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:445
G4double ComputeGroundStateMass(G4int Z, G4int A) const
Definition: G4Fragment.hh:248
void SetSpin(G4double value)
Definition: G4Fragment.hh:414
G4int GetA_asInt() const
Definition: G4Fragment.hh:258
void SetPolarizationFlag(G4bool val)
void SetTwoJMAX(G4int val)
void SetVerbose(G4int val)
virtual G4Fragment * SampleTransition(G4Fragment *nucleus, G4double newExcEnergy, G4double mpRatio, G4int JP1, G4int JP2, G4int MP, G4int shell, G4bool isDiscrete, G4bool isGamma)
G4double LifeTime(size_t i) const
const G4NucLevel * GetLevel(size_t i) const
size_t NearestLevelIndex(G4double energy, size_t index=0) const
G4double NearestLevelEnergy(G4double energy, size_t index=0) const
G4int FloatingLevel(size_t i) const
G4int SpinTwo(size_t i) const
size_t NumberOfTransitions() const
G4double LevelEnergy(size_t i) const
G4int TransitionType(size_t idx) const
Definition: G4NucLevel.hh:121
G4float MultipolarityRatio(size_t idx) const
Definition: G4NucLevel.hh:150
G4float GammaProbability(size_t idx) const
Definition: G4NucLevel.hh:134
size_t NumberOfTransitions() const
Definition: G4NucLevel.hh:108
size_t FinalExcitationIndex(size_t idx) const
Definition: G4NucLevel.hh:113
size_t SampleGammaTransition(G4double rndm) const
Definition: G4NucLevel.hh:158
G4int SampleShell(size_t idx, G4double rndm) const
Definition: G4NucLevel.hh:168
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4NuclearPolarizationStore * GetInstance()
void RemoveMe(G4NuclearPolarization *ptr)
G4NuclearPolarization * FindOrBuild(G4int Z, G4int A, G4double Eexc)
void SetExcitationEnergy(G4double val)
virtual G4double GetFinalLevelEnergy(G4int Z, G4int A, G4double energy) final
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus) final
virtual G4double GetEmissionProbability(G4Fragment *theNucleus) final
virtual void Initialise() final
virtual void SetICM(G4bool)
virtual G4double GetUpperLevelEnergy(G4int Z, G4int A) final
void SetGammaTransition(G4GammaTransition *)
G4PhotonEvaporation(G4GammaTransition *ptr=nullptr)
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)
virtual void RDMForced(G4bool)
virtual G4bool BreakUpChain(G4FragmentVector *theResult, G4Fragment *theNucleus) final
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:225
#define DBL_MAX
Definition: templates.hh:62