Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
34#include "G4ReactionProduct.hh"
35#include "G4Nucleus.hh"
36#include "G4Proton.hh"
37#include "G4Deuteron.hh"
38#include "G4Triton.hh"
39#include "G4Alpha.hh"
40#include "G4ThreeVector.hh"
41#include "G4Poisson.hh"
42#include "G4LorentzVector.hh"
44
46{
47 G4String tString = "/FS/";
48 G4bool dbool;
49 G4ParticleHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
50 G4String filename = aFile.GetName();
51 SetAZMs( A, Z, M, aFile );
52 if(!dbool)
53 {
54 hasAnyData = false;
55 hasFSData = false;
56 hasXsec = false;
57 return;
58 }
59
60 std::istringstream theData(std::ios::in);
62
63 G4int infoType, dataType;
64 hasFSData = false;
65 while (theData >> infoType) // Loop checking, 11.05.2015, T. Koi
66 {
67 hasFSData = true;
68 theData >> dataType;
69 switch(infoType)
70 {
71 case 1:
72 if(dataType==4) theNeutronAngularDis.Init(theData);
73 if(dataType==5) thePromptNeutronEnDis.Init(theData);
74 if(dataType==12) theFinalStatePhotons.InitMean(theData);
75 if(dataType==14) theFinalStatePhotons.InitAngular(theData);
76 if(dataType==15) theFinalStatePhotons.InitEnergies(theData);
77 break;
78 case 2:
79 if(dataType==1) theFinalStateNeutrons.InitMean(theData);
80 break;
81 case 3:
82 if(dataType==1) theFinalStateNeutrons.InitDelayed(theData);
83 if(dataType==5) theDelayedNeutronEnDis.Init(theData);
84 break;
85 case 4:
86 if(dataType==1) theFinalStateNeutrons.InitPrompt(theData);
87 break;
88 case 5:
89 if(dataType==1) theEnergyRelease.Init(theData);
90 break;
91 default:
92 G4cout << "G4ParticleHPFSFissionFS::Init: unknown data type"<<dataType<<G4endl;
93 throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPFSFissionFS::Init: unknown data type");
94 break;
95 }
96 }
97}
98
101 G4double * theDecayConst)
102{
103 G4int i;
105 G4ReactionProduct boosted;
106 boosted.Lorentz( *(fCache.Get().theNeutronRP) , *(fCache.Get().theTarget) );
107 G4double eKinetic = boosted.GetKineticEnergy();
108
109 // Build neutrons
110 G4ReactionProduct * theNeutrons = new G4ReactionProduct[nPrompt+nDelayed];
111 for(i=0; i<nPrompt+nDelayed; ++i)
112 {
113 theNeutrons[i].SetDefinition(G4Neutron::Neutron());
114 }
115
116 // sample energies
117 G4int it, dummy;
118 G4double tempE;
119 for(i=0; i<nPrompt; ++i)
120 {
121 tempE = thePromptNeutronEnDis.Sample(eKinetic, dummy); // energy distribution (file5) always in lab
122 theNeutrons[i].SetKineticEnergy(tempE);
123 }
124 for(i=nPrompt; i<nPrompt+nDelayed; ++i)
125 {
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 {
134 theNeutronAngularDis.SampleAndUpdate(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 {
141 dp->SetDefinition(theNeutrons[i].GetDefinition());
142 dp->SetMomentum(theNeutrons[i].GetMomentum());
143 aResult->push_back(dp);
144 }
145 delete [] theNeutrons;
146
147 return aResult;
148}
149
151{
152 G4double promptNeutronMulti = 0;
153 promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic);
154 G4double delayedNeutronMulti = 0;
155 delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic);
156
157 if(delayedNeutronMulti==0&&promptNeutronMulti==0)
158 {
159 Prompt = 0;
160 delayed = 0;
161 G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
162 all = (G4int)G4Poisson(totalNeutronMulti-off);
163 all += off;
164 }
165 else
166 {
167 Prompt = (G4int)G4Poisson(promptNeutronMulti-off);
168 Prompt += off;
169 delayed = (G4int)G4Poisson(delayedNeutronMulti);
170 all = Prompt+delayed;
171 }
172}
173
175{
176 // sample photons
178 G4ReactionProduct boosted;
179
180 // the photon distributions are in the Nucleus rest frame.
181 boosted.Lorentz( *(fCache.Get().theNeutronRP) , *(fCache.Get().theTarget) );
182 G4double anEnergy = boosted.GetKineticEnergy();
183 temp = theFinalStatePhotons.GetPhotons(anEnergy);
184 if(temp == 0) { return 0; }
185
186 // lorentz transform, and add photons to final state
187 unsigned int i;
189 for(i=0; i<temp->size(); ++i)
190 {
191 // back to lab
192 temp->operator[](i)->Lorentz(*(temp->operator[](i)), -1.* (*(fCache.Get().theTarget)) );
194 theOne->SetDefinition(temp->operator[](i)->GetDefinition());
195 theOne->SetMomentum(temp->operator[](i)->GetMomentum());
196 result->push_back(theOne);
197 delete temp->operator[](i);
198 }
199 delete temp;
200 return result;
201}
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 G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
value_type & Get() const
Definition: G4Cache.hh:315
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
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 *)
G4DynamicParticleVector * ApplyYourself(G4int Prompt, G4int delayed, G4double *decayconst)
void SampleNeutronMult(G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
G4DynamicParticleVector * GetPhotons()
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
void Init(std::istream &aDataFile)
static G4ParticleHPManager * GetInstance()
void GetDataStream(G4String, std::istringstream &iss)
G4ParticleHPDataUsed GetName(G4int A, G4int Z, G4String base, 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)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)