Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPFSFissionFS.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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// P. Arce, June-2014 Conversion neutron_hp to particle_hp
31//
33
34#include "G4Alpha.hh"
35#include "G4Deuteron.hh"
36#include "G4LorentzVector.hh"
37#include "G4Nucleus.hh"
40#include "G4Poisson.hh"
41#include "G4Proton.hh"
42#include "G4ReactionProduct.hh"
43#include "G4ThreeVector.hh"
44#include "G4Triton.hh"
45
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}
100
102 G4double* theDecayConst)
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}
146
148 G4double eKinetic, G4int off)
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}
169
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< G4DynamicParticle * > G4DynamicParticleVector
#define M(row, col)
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
std::vector< G4ReactionProduct * > G4ReactionProductVector
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
value_type & Get() const
Definition G4Cache.hh:315
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
void Init(std::istream &aDataFile)
void SampleAndUpdate(G4ReactionProduct &anIncidentParticle)
G4double Sample(G4double anEnergy, G4int &it)
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *) override
G4DynamicParticleVector * ApplyYourself(G4int Prompt, G4int delayed, G4double *decayconst)
void SampleNeutronMult(G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
G4DynamicParticleVector * GetPhotons()
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)
G4double GetMean(G4double anEnergy)
void InitMean(std::istream &aDataFile)
void InitDelayed(std::istream &aDataFile)
G4double GetPrompt(G4double anEnergy)
void InitPrompt(std::istream &aDataFile)
G4double GetDelayed(G4double anEnergy)
void InitEnergies(std::istream &aDataFile)
G4ReactionProductVector * GetPhotons(G4double anEnergy)
void InitAngular(std::istream &aDataFile)
G4bool InitMean(std::istream &aDataFile)
G4double GetKineticEnergy() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)