Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPFSFissionFS.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 //
31 #include "G4NeutronHPManager.hh"
32 #include "G4ReactionProduct.hh"
33 #include "G4Nucleus.hh"
34 #include "G4Proton.hh"
35 #include "G4Deuteron.hh"
36 #include "G4Triton.hh"
37 #include "G4Alpha.hh"
38 #include "G4ThreeVector.hh"
39 #include "G4Poisson.hh"
40 #include "G4LorentzVector.hh"
41 #include "G4NeutronHPDataUsed.hh"
42 
44  {
45  G4String tString = "/FS/";
46  G4bool dbool;
47  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
48  G4String filename = aFile.GetName();
49  SetAZMs( A, Z, M, aFile );
50  if(!dbool)
51  {
52  hasAnyData = false;
53  hasFSData = false;
54  hasXsec = false;
55  return;
56  }
57  //std::ifstream theData(filename, std::ios::in);
58  std::istringstream theData(std::ios::in);
60 
61  // here it comes
62  G4int infoType, dataType;
63  hasFSData = false;
64  while (theData >> infoType)
65  {
66  hasFSData = true;
67  theData >> dataType;
68  switch(infoType)
69  {
70  case 1:
71  if(dataType==4) theNeutronAngularDis.Init(theData);
72  if(dataType==5) thePromptNeutronEnDis.Init(theData);
73  if(dataType==12) theFinalStatePhotons.InitMean(theData);
74  if(dataType==14) theFinalStatePhotons.InitAngular(theData);
75  if(dataType==15) theFinalStatePhotons.InitEnergies(theData);
76  break;
77  case 2:
78  if(dataType==1) theFinalStateNeutrons.InitMean(theData);
79  break;
80  case 3:
81  if(dataType==1) theFinalStateNeutrons.InitDelayed(theData);
82  if(dataType==5) theDelayedNeutronEnDis.Init(theData);
83  break;
84  case 4:
85  if(dataType==1) theFinalStateNeutrons.InitPrompt(theData);
86  break;
87  case 5:
88  if(dataType==1) theEnergyRelease.Init(theData);
89  break;
90  default:
91  G4cout << "G4NeutronHPFSFissionFS::Init: unknown data type"<<dataType<<G4endl;
92  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPFSFissionFS::Init: unknown data type");
93  break;
94  }
95  }
96  targetMass = theFinalStateNeutrons.GetTargetMass();
97  //theData.close();
98  }
99 
100 
102  G4int nDelayed, G4double * theDecayConst)
103  {
104  G4int i;
106  G4ReactionProduct boosted;
107  boosted.Lorentz(theNeutron, theTarget);
108  G4double eKinetic = boosted.GetKineticEnergy();
109 
110 // Build neutrons
111  G4ReactionProduct * theNeutrons = new G4ReactionProduct[nPrompt+nDelayed];
112  for(i=0; i<nPrompt+nDelayed; i++)
113  {
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  {
122  tempE = 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  {
127  theNeutrons[i].SetKineticEnergy(theDelayedNeutronEnDis.Sample(eKinetic, it)); // dito
128  if(it==0) theNeutrons[i].SetKineticEnergy(thePromptNeutronEnDis.Sample(eKinetic, dummy));
129  theDecayConst[i-nPrompt] = theFinalStateNeutrons.GetDecayConstant(it); // this is returned
130  }
131 
132 // sample neutron angular distribution
133  for(i=0; i<nPrompt+nDelayed; i++)
134  {
135  theNeutronAngularDis.SampleAndUpdate(theNeutrons[i]); // angular comes back in lab automatically
136  }
137 
138 // already in lab. Add neutrons to dynamic particle vector
139  for(i=0; i<nPrompt+nDelayed; i++)
140  {
142  dp->SetDefinition(theNeutrons[i].GetDefinition());
143  dp->SetMomentum(theNeutrons[i].GetMomentum());
144  aResult->push_back(dp);
145  }
146  delete [] theNeutrons;
147 // return the result
148  return aResult;
149  }
150 
151 void G4NeutronHPFSFissionFS::SampleNeutronMult(G4int&all, G4int&Prompt, G4int&delayed, G4double eKinetic, G4int off)
152 {
153  G4double promptNeutronMulti = 0;
154  promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic);
155  G4double delayedNeutronMulti = 0;
156  delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic);
157 
158  if(delayedNeutronMulti==0&&promptNeutronMulti==0)
159  {
160  Prompt = 0;
161  delayed = 0;
162  G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
163  all = G4Poisson(totalNeutronMulti-off);
164  all += off;
165  }
166  else
167  {
168  Prompt = G4Poisson(promptNeutronMulti-off);
169  Prompt += off;
170  delayed = G4Poisson(delayedNeutronMulti);
171  all = Prompt+delayed;
172  }
173 }
174 
176 {
177 // sample photons
179  G4ReactionProduct boosted;
180 // the photon distributions are in the Nucleus rest frame.
181  boosted.Lorentz(theNeutron, 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.*theTarget);
193  G4DynamicParticle * theOne = new G4DynamicParticle;
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 }
G4ReactionProductVector * GetPhotons(G4double anEnergy)
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
void SetMomentum(const G4ThreeVector &momentum)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void Init(std::istream &aDataFile)
void SetKineticEnergy(const G4double en)
static G4NeutronHPManager * GetInstance()
void InitDelayed(std::istream &aDataFile)
G4double Sample(G4double anEnergy, G4int &it)
void GetDataStream(G4String, std::istringstream &iss)
G4NeutronHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
int G4int
Definition: G4Types.hh:78
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
G4double GetMean(G4double anEnergy)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void InitAngular(std::istream &aDataFile)
G4GLOB_DLL std::ostream G4cout
G4DynamicParticleVector * ApplyYourself(G4int Prompt, G4int delayed, G4double *decayconst)
bool G4bool
Definition: G4Types.hh:79
void InitEnergies(std::istream &aDataFile)
G4DynamicParticleVector * GetPhotons()
std::vector< G4DynamicParticle * > G4DynamicParticleVector
void Init(std::istream &aDataFile)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetKineticEnergy() const
void SampleAndUpdate(G4ReactionProduct &aNeutron)
G4double GetDelayed(G4double anEnergy)
G4bool InitMean(std::istream &aDataFile)
void SampleNeutronMult(G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
void InitPrompt(std::istream &aDataFile)
#define G4endl
Definition: G4ios.hh:61
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4NeutronHPDataUsed used)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetPrompt(G4double anEnergy)
void InitMean(std::istream &aDataFile)