Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPFissionFS Class Reference

#include <G4ParticleHPFissionFS.hh>

+ Inheritance diagram for G4ParticleHPFissionFS:

Public Member Functions

 G4ParticleHPFissionFS ()
 
 ~G4ParticleHPFissionFS () override=default
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *) override
 
G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack) override
 
G4ParticleHPFinalStateNew () override
 
- Public Member Functions inherited from G4ParticleHPFinalState
 G4ParticleHPFinalState ()
 
virtual ~G4ParticleHPFinalState ()
 
void Init (G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *p)
 
G4bool HasXsec ()
 
G4bool HasFSData ()
 
G4bool HasAnyData ()
 
virtual G4double GetXsec (G4double)
 
virtual G4ParticleHPVectorGetXsec ()
 
void SetA_Z (G4double anA, G4double aZ, G4int aM=0)
 
G4double GetZ ()
 
G4double GetN ()
 
G4double GetA ()
 
G4int GetM ()
 
void SetAZMs (G4ParticleHPDataUsed used)
 
void SetAZMs (G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
 
void SetProjectile (G4ParticleDefinition *projectile)
 
G4ParticleHPFinalStateoperator= (const G4ParticleHPFinalState &right)=delete
 
 G4ParticleHPFinalState (const G4ParticleHPFinalState &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4ParticleHPFinalState
void adjust_final_state (G4LorentzVector)
 
- Protected Attributes inherited from G4ParticleHPFinalState
G4ParticleDefinitiontheProjectile {nullptr}
 
G4ParticleHPManagerfManager
 
G4IonTableionTable
 
G4int theBaseA {0}
 
G4int theBaseZ {0}
 
G4int theBaseM {0}
 
G4int theNDLDataZ {0}
 
G4int theNDLDataA {0}
 
G4int theNDLDataM {0}
 
G4int secID {-1}
 
G4bool hasXsec {true}
 
G4bool hasFSData {true}
 
G4bool hasAnyData {true}
 
G4ParticleHPNames theNames
 
G4Cache< G4HadFinalState * > theResult
 

Detailed Description

Definition at line 44 of file G4ParticleHPFissionFS.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFissionFS()

G4ParticleHPFissionFS::G4ParticleHPFissionFS ( )

Definition at line 45 of file G4ParticleHPFissionFS.cc.

46{
47 secID = G4PhysicsModelCatalog::GetModelID("model_NeutronHPFission");
48 hasXsec = false;
49 produceFissionFragments = false;
50}
static G4int GetModelID(const G4int modelIndex)

Referenced by New().

◆ ~G4ParticleHPFissionFS()

G4ParticleHPFissionFS::~G4ParticleHPFissionFS ( )
overridedefault

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ParticleHPFissionFS::ApplyYourself ( const G4HadProjectile & theTrack)
overridevirtual

Reimplemented from G4ParticleHPFinalState.

Definition at line 72 of file G4ParticleHPFissionFS.cc.

73{
74 // Because it may change by UI command
76
77 // prepare neutron
78 if (theResult.Get() == nullptr) theResult.Put(new G4HadFinalState);
79 theResult.Get()->Clear();
80 G4double eKinetic = theTrack.GetKineticEnergy();
81 const G4HadProjectile* incidentParticle = &theTrack;
82 G4ReactionProduct theNeutron(
83 const_cast<G4ParticleDefinition*>(incidentParticle->GetDefinition()));
84 theNeutron.SetMomentum(incidentParticle->Get4Momentum().vect());
85 theNeutron.SetKineticEnergy(eKinetic);
86
87 // prepare target
88 G4Nucleus aNucleus;
89 G4ReactionProduct theTarget;
90 G4double targetMass = theFS.GetMass();
91 G4ThreeVector neuVelo =
92 (1. / incidentParticle->GetDefinition()->GetPDGMass()) * theNeutron.GetMomentum();
93 theTarget =
94 aNucleus.GetBiasedThermalNucleus(targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
95 theTarget.SetDefinition(
96 G4IonTable::GetIonTable()->GetIon(G4int(theBaseZ), G4int(theBaseA), 0.0)); // TESTPHP
97
98 // set neutron and target in the FS classes
99 theFS.SetNeutronRP(theNeutron);
100 theFS.SetTarget(theTarget);
101 theFC.SetNeutronRP(theNeutron);
102 theFC.SetTarget(theTarget);
103 theSC.SetNeutronRP(theNeutron);
104 theSC.SetTarget(theTarget);
105 theTC.SetNeutronRP(theNeutron);
106 theTC.SetTarget(theTarget);
107 theLC.SetNeutronRP(theNeutron);
108 theLC.SetTarget(theTarget);
109 theFF.SetNeutronRP(theNeutron);
110 theFF.SetTarget(theTarget);
111
112 // boost to target rest system and decide on channel.
113 theNeutron.Lorentz(theNeutron, -1 * theTarget);
114
115 // dice the photons
116
117 G4DynamicParticleVector* thePhotons;
118 thePhotons = theFS.GetPhotons();
119
120 // select the FS in charge
121
122 eKinetic = theNeutron.GetKineticEnergy();
123 G4double xSec[4];
124 xSec[0] = theFC.GetXsec(eKinetic);
125 xSec[1] = xSec[0] + theSC.GetXsec(eKinetic);
126 xSec[2] = xSec[1] + theTC.GetXsec(eKinetic);
127 xSec[3] = xSec[2] + theLC.GetXsec(eKinetic);
128 G4int it;
129 unsigned int i = 0;
130 G4double random = G4UniformRand();
131 if (xSec[3] == 0) {
132 it = -1;
133 }
134 else {
135 for (i = 0; i < 4; i++) {
136 it = i;
137 if (random < xSec[i] / xSec[3]) break;
138 }
139 }
140
141 // dice neutron multiplicities, energies and momenta in Lab. @@
142 // no energy conservation on an event-to-event basis. we rely on the data to be ok. @@
143 // also for mean, we rely on the consistancy of the data. @@
144
145 G4int Prompt = 0, delayed = 0, all = 0;
146 G4DynamicParticleVector* theNeutrons = nullptr;
147 switch (it) // check logic, and ask, if partials can be assumed to correspond to individual
148 // particles @@@
149 {
150 case 0:
151 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0);
152 if (Prompt == 0 && delayed == 0) Prompt = all;
153 theNeutrons = theFC.ApplyYourself(Prompt); // delayed always in FS
154 // take 'U' into account explicitly (see 5.4) in the sampling of energy @@@@
155 break;
156 case 1:
157 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 1);
158 if (Prompt == 0 && delayed == 0) Prompt = all;
159 theNeutrons = theSC.ApplyYourself(Prompt); // delayed always in FS, off done in FSFissionFS
160 break;
161 case 2:
162 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 2);
163 if (Prompt == 0 && delayed == 0) Prompt = all;
164 theNeutrons = theTC.ApplyYourself(Prompt); // delayed always in FS
165 break;
166 case 3:
167 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 3);
168 if (Prompt == 0 && delayed == 0) Prompt = all;
169 theNeutrons = theLC.ApplyYourself(Prompt); // delayed always in FS
170 break;
171 default:
172 break;
173 }
174
175 // dice delayed neutrons and photons, and fallback
176 // for Prompt in case channel had no FS data; add all paricles to FS.
177
178 if (produceFissionFragments) delayed = 0;
179
180 G4double* theDecayConstants;
181
182 if (theNeutrons != nullptr) {
183 theDecayConstants = new G4double[delayed];
184 for (i = 0; i < theNeutrons->size(); ++i) {
185 theResult.Get()->AddSecondary(theNeutrons->operator[](i), secID);
186 }
187 delete theNeutrons;
188
189 G4DynamicParticleVector* theDelayed = nullptr;
190 theDelayed = theFS.ApplyYourself(0, delayed, theDecayConstants);
191 for (i = 0; i < theDelayed->size(); i++) {
192 G4double time = -G4Log(G4UniformRand()) / theDecayConstants[i];
193 time += theTrack.GetGlobalTime();
194 theResult.Get()->AddSecondary(theDelayed->operator[](i), secID);
196 }
197 delete theDelayed;
198 }
199 else {
200 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0);
201 theDecayConstants = new G4double[delayed];
202 if (Prompt == 0 && delayed == 0) Prompt = all;
203 theNeutrons = theFS.ApplyYourself(Prompt, delayed, theDecayConstants);
204 G4int i0;
205 for (i0 = 0; i0 < Prompt; ++i0) {
206 theResult.Get()->AddSecondary(theNeutrons->operator[](i0), secID);
207 }
208
209 for (i0 = Prompt; i0 < Prompt + delayed; ++i0) {
210 // Protect against the very rare case of division by zero
211 G4double time = 0.0;
212 if (theDecayConstants[i0 - Prompt] > 1.0e-30) {
213 time = -G4Log(G4UniformRand()) / theDecayConstants[i0 - Prompt];
214 }
215 else {
217 ed << " theDecayConstants[i0-Prompt]=" << theDecayConstants[i0 - Prompt]
218 << " -> cannot sample the time : set it to 0.0 !" << G4endl;
219 G4Exception("G4ParticleHPFissionFS::ApplyYourself ", "HAD_FISSIONHP_001", JustWarning, ed);
220 }
221
222 time += theTrack.GetGlobalTime();
223 theResult.Get()->AddSecondary(theNeutrons->operator[](i0), secID);
225 }
226 delete theNeutrons;
227 }
228 delete[] theDecayConstants;
229
230 std::size_t nPhotons = 0;
231 if (thePhotons != nullptr) {
232 nPhotons = thePhotons->size();
233 for (i = 0; i < nPhotons; ++i) {
234 theResult.Get()->AddSecondary(thePhotons->operator[](i), secID);
235 }
236 delete thePhotons;
237 }
238
239 // finally deal with local energy depositions.
240
241 G4ParticleHPFissionERelease* theERelease = theFS.GetEnergyRelease();
242 G4double eDepByFragments = theERelease->GetFragmentKinetic();
243 // theResult.SetLocalEnergyDeposit(eDepByFragments);
244 if (!produceFissionFragments) theResult.Get()->SetLocalEnergyDeposit(eDepByFragments);
245 // clean up the primary neutron
247
248 if (produceFissionFragments) {
249 G4int fragA_Z = 0;
250 G4int fragA_A = 0;
251 G4int fragA_M = 0;
252 // System is traget rest!
253 theFF.GetAFissionFragment(eKinetic, fragA_Z, fragA_A, fragA_M);
254 if (0 == fragA_A) { return theResult.Get(); }
255 G4int fragB_Z = (G4int)theBaseZ - fragA_Z;
256 G4int fragB_A = (G4int)theBaseA - fragA_A - Prompt;
257
259 // Excitation energy is not taken into account
260 G4ParticleDefinition* pdA = pt->GetIon(fragA_Z, fragA_A, 0.0);
261 G4ParticleDefinition* pdB = pt->GetIon(fragB_Z, fragB_A, 0.0);
262
263 // Isotropic Distribution
264 G4double phi = twopi * G4UniformRand();
265 // Bug #1745 DHW G4double theta = pi*G4UniformRand();
266 G4double costheta = 2. * G4UniformRand() - 1.;
267 G4double theta = std::acos(costheta);
268 G4double sinth = std::sin(theta);
269 G4ThreeVector direction(sinth * std::cos(phi), sinth * std::sin(phi), costheta);
270
271 // Just use ENDF value for this
272 G4double ER = eDepByFragments;
273 G4double ma = pdA->GetPDGMass();
274 G4double mb = pdB->GetPDGMass();
275 G4double EA = ER / (1 + ma / mb);
276 G4double EB = ER - EA;
277 auto dpA = new G4DynamicParticle(pdA, direction, EA);
278 auto dpB = new G4DynamicParticle(pdB, -direction, EB);
281 }
282
283 return theResult.Get();
284}
std::vector< G4DynamicParticle * > G4DynamicParticleVector
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ stopAndKill
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector vect() const
value_type & Get() const
Definition G4Cache.hh:315
void Put(const value_type &val) const
Definition G4Cache.hh:321
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void SetLocalEnergyDeposit(G4double aE)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4IonTable * GetIonTable()
G4double GetTemperature() const
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition G4Nucleus.cc:119
G4DynamicParticleVector * ApplyYourself(G4int nNeutrons)
void GetAFissionFragment(G4double, G4int &, G4int &, G4int &)
G4DynamicParticleVector * ApplyYourself(G4int Prompt, G4int delayed, G4double *decayconst)
void SampleNeutronMult(G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
G4ParticleHPFissionERelease * GetEnergyRelease()
void SetNeutronRP(const G4ReactionProduct &aNeutron)
void SetTarget(const G4ReactionProduct &aTarget)
G4DynamicParticleVector * GetPhotons()
G4Cache< G4HadFinalState * > theResult
G4double GetXsec(G4double anEnergy) override
void SetTarget(const G4ReactionProduct &aTarget)
void SetNeutronRP(const G4ReactionProduct &aNeutron)
G4DynamicParticleVector * ApplyYourself(G4int NNeutrons)
G4bool GetProduceFissionFragments() const
static G4ParticleHPManager * GetInstance()
G4DynamicParticleVector * ApplyYourself(G4int NNeutrons)
G4DynamicParticleVector * ApplyYourself(G4int NNeutrons)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)

◆ Init()

void G4ParticleHPFissionFS::Init ( G4double A,
G4double Z,
G4int M,
G4String & dirName,
G4String & aFSType,
G4ParticleDefinition * projectile )
overridevirtual

Implements G4ParticleHPFinalState.

Definition at line 52 of file G4ParticleHPFissionFS.cc.

54{
55 theFS.Init(A, Z, M, dirName, aFSType, projectile);
56 theFC.Init(A, Z, M, dirName, aFSType, projectile);
57 theSC.Init(A, Z, M, dirName, aFSType, projectile);
58 theTC.Init(A, Z, M, dirName, aFSType, projectile);
59 theLC.Init(A, Z, M, dirName, aFSType, projectile);
60
61 theFF.Init(A, Z, M, dirName, aFSType, projectile);
62 if (G4ParticleHPManager::GetInstance()->GetProduceFissionFragments() && theFF.HasFSData()) {
63 G4cout << "Fission fragment production is now activated in HP package for "
64 << "Z = " << (G4int)Z << ", A = " << (G4int)A << G4endl;
65 G4cout << "As currently modeled this option precludes production of delayed neutrons from "
66 "fission fragments."
67 << G4endl;
68 produceFissionFragments = true;
69 }
70}
#define M(row, col)
const G4double A[17]
G4GLOB_DLL std::ostream G4cout
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile) override
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *) override
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *) override
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile) override
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile) override
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile) override

◆ New()

G4ParticleHPFinalState * G4ParticleHPFissionFS::New ( )
inlineoverridevirtual

Implements G4ParticleHPFinalState.

Definition at line 52 of file G4ParticleHPFissionFS.hh.

53 {
54 auto theNew = new G4ParticleHPFissionFS;
55 return theNew;
56 }

The documentation for this class was generated from the following files: