Geant4-11
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 //std::ifstream theData(filename, std::ios::in);
60 std::istringstream theData(std::ios::in);
62
63 // here it comes
64 G4int infoType, dataType;
65 hasFSData = false;
66 while (theData >> infoType) // Loop checking, 11.05.2015, T. Koi
67 {
68 hasFSData = true;
69 theData >> dataType;
70 switch(infoType)
71 {
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__, "G4ParticleHPFSFissionFS::Init: unknown data type");
95 break;
96 }
97 }
98 //targetMass = theFinalStateNeutrons.GetTargetMass();
99 //theData.close();
100 }
101
102
104 G4int nDelayed, G4double * theDecayConst)
105 {
106 G4int i;
108 G4ReactionProduct boosted;
109 boosted.Lorentz( *(fCache.Get().theNeutronRP) , *(fCache.Get().theTarget) );
110 G4double eKinetic = boosted.GetKineticEnergy();
111
112// Build neutrons
113 G4ReactionProduct * theNeutrons = new G4ReactionProduct[nPrompt+nDelayed];
114 for(i=0; i<nPrompt+nDelayed; i++)
115 {
116 theNeutrons[i].SetDefinition(G4Neutron::Neutron());
117 }
118
119// sample energies
120 G4int it, dummy;
121 G4double tempE;
122 for(i=0; i<nPrompt; i++)
123 {
124 tempE = thePromptNeutronEnDis.Sample(eKinetic, dummy); // energy distribution (file5) always in lab
125 theNeutrons[i].SetKineticEnergy(tempE);
126 }
127 for(i=nPrompt; i<nPrompt+nDelayed; i++)
128 {
129 theNeutrons[i].SetKineticEnergy(theDelayedNeutronEnDis.Sample(eKinetic, it)); // dito
130 if(it==0) theNeutrons[i].SetKineticEnergy(thePromptNeutronEnDis.Sample(eKinetic, dummy));
131 theDecayConst[i-nPrompt] = theFinalStateNeutrons.GetDecayConstant(it); // this is returned
132 }
133
134// sample neutron angular distribution
135 for(i=0; i<nPrompt+nDelayed; i++)
136 {
137 theNeutronAngularDis.SampleAndUpdate(theNeutrons[i]); // angular comes back in lab automatically
138 }
139
140// already in lab. Add neutrons to dynamic particle vector
141 for(i=0; i<nPrompt+nDelayed; i++)
142 {
144 dp->SetDefinition(theNeutrons[i].GetDefinition());
145 dp->SetMomentum(theNeutrons[i].GetMomentum());
146 aResult->push_back(dp);
147 }
148 delete [] theNeutrons;
149// return the result
150 return aResult;
151 }
152
154{
155 G4double promptNeutronMulti = 0;
156 promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic);
157 G4double delayedNeutronMulti = 0;
158 delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic);
159
160 if(delayedNeutronMulti==0&&promptNeutronMulti==0)
161 {
162 Prompt = 0;
163 delayed = 0;
164 G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
165 all = G4Poisson(totalNeutronMulti-off);
166 all += off;
167 }
168 else
169 {
170 Prompt = G4Poisson(promptNeutronMulti-off);
171 Prompt += off;
172 delayed = G4Poisson(delayedNeutronMulti);
173 all = Prompt+delayed;
174 }
175}
176
178{
179// sample photons
181 G4ReactionProduct boosted;
182// the photon distributions are in the Nucleus rest frame.
183 boosted.Lorentz( *(fCache.Get().theNeutronRP) , *(fCache.Get().theTarget) );
184 G4double anEnergy = boosted.GetKineticEnergy();
185 temp = theFinalStatePhotons.GetPhotons(anEnergy);
186 if(temp == 0) { return 0; }
187
188// lorentz transform, and add photons to final state
189 unsigned int i;
191 for(i=0; i<temp->size(); i++)
192 {
193 // back to lab
194 temp->operator[](i)->Lorentz(*(temp->operator[](i)), -1.* (*(fCache.Get().theTarget)) );
196 theOne->SetDefinition(temp->operator[](i)->GetDefinition());
197 theOne->SetMomentum(temp->operator[](i)->GetMomentum());
198 result->push_back(theOne);
199 delete temp->operator[](i);
200 }
201 delete temp;
202 return result;
203}
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
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()
G4ParticleHPAngular theNeutronAngularDis
G4ParticleHPFissionERelease theEnergyRelease
G4ParticleHPPhotonDist theFinalStatePhotons
G4ParticleHPParticleYield theFinalStateNeutrons
G4ParticleHPEnergyDistribution thePromptNeutronEnDis
G4ParticleHPEnergyDistribution theDelayedNeutronEnDis
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)