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

#include <G4ParticleHPFSFissionFS.hh>

+ Inheritance diagram for G4ParticleHPFSFissionFS:

Public Member Functions

 G4ParticleHPFSFissionFS ()
 
 ~G4ParticleHPFSFissionFS () override=default
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *) override
 
G4DynamicParticleVectorApplyYourself (G4int Prompt, G4int delayed, G4double *decayconst)
 
G4ParticleHPFinalStateNew () override
 
G4double GetMass ()
 
void SampleNeutronMult (G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
 
void SetNeutronRP (const G4ReactionProduct &aNeutron)
 
void SetTarget (const G4ReactionProduct &aTarget)
 
G4DynamicParticleVectorGetPhotons ()
 
G4ParticleHPFissionEReleaseGetEnergyRelease ()
 
- 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 45 of file G4ParticleHPFSFissionFS.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFSFissionFS()

G4ParticleHPFSFissionFS::G4ParticleHPFSFissionFS ( )
inline

Definition at line 55 of file G4ParticleHPFSFissionFS.hh.

Referenced by New().

◆ ~G4ParticleHPFSFissionFS()

G4ParticleHPFSFissionFS::~G4ParticleHPFSFissionFS ( )
overridedefault

Member Function Documentation

◆ ApplyYourself()

G4DynamicParticleVector * G4ParticleHPFSFissionFS::ApplyYourself ( G4int Prompt,
G4int delayed,
G4double * decayconst )

Definition at line 101 of file G4ParticleHPFSFissionFS.cc.

103{
104 G4int i;
105 auto aResult = new G4DynamicParticleVector;
106 G4ReactionProduct boosted;
107 boosted.Lorentz(*(fCache.Get().theNeutronRP), *(fCache.Get().theTarget));
108 G4double eKinetic = boosted.GetKineticEnergy();
109
110 // Build neutrons
111 std::vector<G4ReactionProduct> theNeutrons;
112 for (i = 0; i < nPrompt + nDelayed; ++i) {
113 theNeutrons.emplace_back();
114 theNeutrons[i].SetDefinition(G4Neutron::Neutron());
115 }
116
117 // sample energies
118 G4int it, dummy;
119 G4double tempE;
120 for (i = 0; i < nPrompt; ++i) {
121 tempE =
122 thePromptNeutronEnDis.Sample(eKinetic, dummy); // energy distribution (file5) always in lab
123 theNeutrons[i].SetKineticEnergy(tempE);
124 }
125 for (i = nPrompt; i < nPrompt + nDelayed; ++i) {
126 theNeutrons[i].SetKineticEnergy(theDelayedNeutronEnDis.Sample(eKinetic, it)); // dito
127 if (it == 0) theNeutrons[i].SetKineticEnergy(thePromptNeutronEnDis.Sample(eKinetic, dummy));
128 theDecayConst[i - nPrompt] = theFinalStateNeutrons.GetDecayConstant(it); // this is returned
129 }
130
131 // sample neutron angular distribution
132 for (i = 0; i < nPrompt + nDelayed; ++i) {
133 theNeutronAngularDis.SampleAndUpdate(
134 theNeutrons[i]); // angular comes back in lab automatically
135 }
136
137 // already in lab. Add neutrons to dynamic particle vector
138 for (i = 0; i < nPrompt + nDelayed; ++i) {
139 auto dp = new G4DynamicParticle;
140 dp->SetDefinition(theNeutrons[i].GetDefinition());
141 dp->SetMomentum(theNeutrons[i].GetMomentum());
142 aResult->push_back(dp);
143 }
144 return aResult;
145}
std::vector< G4DynamicParticle * > G4DynamicParticleVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
value_type & Get() const
Definition G4Cache.hh:315
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
void SampleAndUpdate(G4ReactionProduct &anIncidentParticle)
G4double Sample(G4double anEnergy, G4int &it)
G4double GetKineticEnergy() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)

Referenced by G4ParticleHPFissionFS::ApplyYourself().

◆ GetEnergyRelease()

G4ParticleHPFissionERelease * G4ParticleHPFSFissionFS::GetEnergyRelease ( )
inline

Definition at line 87 of file G4ParticleHPFSFissionFS.hh.

87{ return &theEnergyRelease; }

Referenced by G4ParticleHPFissionFS::ApplyYourself().

◆ GetMass()

G4double G4ParticleHPFSFissionFS::GetMass ( )
inline

Definition at line 69 of file G4ParticleHPFSFissionFS.hh.

69{ return theFinalStateNeutrons.GetTargetMass(); }

Referenced by G4ParticleHPFissionFS::ApplyYourself().

◆ GetPhotons()

G4DynamicParticleVector * G4ParticleHPFSFissionFS::GetPhotons ( )

Definition at line 170 of file G4ParticleHPFSFissionFS.cc.

171{
172 // sample photons
174 G4ReactionProduct boosted;
175
176 // the photon distributions are in the Nucleus rest frame.
177 boosted.Lorentz(*(fCache.Get().theNeutronRP), *(fCache.Get().theTarget));
178 G4double anEnergy = boosted.GetKineticEnergy();
179 temp = theFinalStatePhotons.GetPhotons(anEnergy);
180 if (temp == nullptr) {
181 return nullptr;
182 }
183
184 // lorentz transform, and add photons to final state
185 unsigned int i;
186 auto result = new G4DynamicParticleVector;
187 for (i = 0; i < temp->size(); ++i) {
188 // back to lab
189 temp->operator[](i)->Lorentz(*(temp->operator[](i)), -1. * (*(fCache.Get().theTarget)));
190 auto theOne = new G4DynamicParticle;
191 theOne->SetDefinition(temp->operator[](i)->GetDefinition());
192 theOne->SetMomentum(temp->operator[](i)->GetMomentum());
193 result->push_back(theOne);
194 delete temp->operator[](i);
195 }
196 delete temp;
197 return result;
198}
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4ReactionProductVector * GetPhotons(G4double anEnergy)

Referenced by G4ParticleHPFissionFS::ApplyYourself().

◆ Init()

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

Implements G4ParticleHPFinalState.

Definition at line 46 of file G4ParticleHPFSFissionFS.cc.

48{
49 G4String tString = "/FS/";
50 G4bool dbool;
52 theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
53 G4String filename = aFile.GetName();
54 SetAZMs(A, Z, M, aFile);
55 if (!dbool) {
56 hasAnyData = false;
57 hasFSData = false;
58 hasXsec = false;
59 return;
60 }
61
62 std::istringstream theData(std::ios::in);
64
65 G4int infoType, dataType;
66 hasFSData = false;
67 while (theData >> infoType) // Loop checking, 11.05.2015, T. Koi
68 {
69 hasFSData = true;
70 theData >> dataType;
71 switch (infoType) {
72 case 1:
73 if (dataType == 4) theNeutronAngularDis.Init(theData);
74 if (dataType == 5) thePromptNeutronEnDis.Init(theData);
75 if (dataType == 12) theFinalStatePhotons.InitMean(theData);
76 if (dataType == 14) theFinalStatePhotons.InitAngular(theData);
77 if (dataType == 15) theFinalStatePhotons.InitEnergies(theData);
78 break;
79 case 2:
80 if (dataType == 1) theFinalStateNeutrons.InitMean(theData);
81 break;
82 case 3:
83 if (dataType == 1) theFinalStateNeutrons.InitDelayed(theData);
84 if (dataType == 5) theDelayedNeutronEnDis.Init(theData);
85 break;
86 case 4:
87 if (dataType == 1) theFinalStateNeutrons.InitPrompt(theData);
88 break;
89 case 5:
90 if (dataType == 1) theEnergyRelease.Init(theData);
91 break;
92 default:
93 G4cout << "G4ParticleHPFSFissionFS::Init: unknown data type" << dataType << G4endl;
94 throw G4HadronicException(__FILE__, __LINE__,
95 "G4ParticleHPFSFissionFS::Init: unknown data type");
96 break;
97 }
98 }
99}
#define M(row, col)
bool G4bool
Definition G4Types.hh:86
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void Init(std::istream &aDataFile)
void SetAZMs(G4ParticleHPDataUsed used)
void Init(std::istream &aDataFile)
void GetDataStream(const G4String &, std::istringstream &iss)
static G4ParticleHPManager * GetInstance()
G4ParticleHPDataUsed GetName(G4int A, G4int Z, const G4String &base, const G4String &rest, G4bool &active)
void InitMean(std::istream &aDataFile)
void InitDelayed(std::istream &aDataFile)
void InitPrompt(std::istream &aDataFile)
void InitEnergies(std::istream &aDataFile)
void InitAngular(std::istream &aDataFile)
G4bool InitMean(std::istream &aDataFile)

Referenced by G4ParticleHPFissionFS::Init().

◆ New()

G4ParticleHPFinalState * G4ParticleHPFSFissionFS::New ( )
inlineoverridevirtual

Implements G4ParticleHPFinalState.

Definition at line 63 of file G4ParticleHPFSFissionFS.hh.

64 {
65 auto theNew = new G4ParticleHPFSFissionFS;
66 return theNew;
67 }

◆ SampleNeutronMult()

void G4ParticleHPFSFissionFS::SampleNeutronMult ( G4int & all,
G4int & Prompt,
G4int & delayed,
G4double energy,
G4int off )

Definition at line 147 of file G4ParticleHPFSFissionFS.cc.

149{
150 G4double promptNeutronMulti = 0;
151 promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic);
152 G4double delayedNeutronMulti = 0;
153 delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic);
154
155 if (delayedNeutronMulti == 0 && promptNeutronMulti == 0) {
156 Prompt = 0;
157 delayed = 0;
158 G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
159 all = (G4int)G4Poisson(totalNeutronMulti - off);
160 all += off;
161 }
162 else {
163 Prompt = (G4int)G4Poisson(promptNeutronMulti - off);
164 Prompt += off;
165 delayed = (G4int)G4Poisson(delayedNeutronMulti);
166 all = Prompt + delayed;
167 }
168}
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
G4double GetMean(G4double anEnergy)
G4double GetPrompt(G4double anEnergy)
G4double GetDelayed(G4double anEnergy)

Referenced by G4ParticleHPFissionFS::ApplyYourself().

◆ SetNeutronRP()

void G4ParticleHPFSFissionFS::SetNeutronRP ( const G4ReactionProduct & aNeutron)
inline

Definition at line 73 of file G4ParticleHPFSFissionFS.hh.

74 {
75 fCache.Get().theNeutronRP = &aNeutron;
76 theNeutronAngularDis.SetProjectileRP(aNeutron);
77 }
void SetProjectileRP(const G4ReactionProduct &anIncidentParticleRP)

Referenced by G4ParticleHPFissionFS::ApplyYourself().

◆ SetTarget()

void G4ParticleHPFSFissionFS::SetTarget ( const G4ReactionProduct & aTarget)
inline

Definition at line 79 of file G4ParticleHPFSFissionFS.hh.

80 {
81 fCache.Get().theTarget = &aTarget;
82 theNeutronAngularDis.SetTarget(aTarget);
83 }
void SetTarget(const G4ReactionProduct &aTarget)

Referenced by G4ParticleHPFissionFS::ApplyYourself().


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