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