Geant4  9.6.p02
 All Classes 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 "G4ReactionProduct.hh"
32 #include "G4Nucleus.hh"
33 #include "G4Proton.hh"
34 #include "G4Deuteron.hh"
35 #include "G4Triton.hh"
36 #include "G4Alpha.hh"
37 #include "G4ThreeVector.hh"
38 #include "G4Poisson.hh"
39 #include "G4LorentzVector.hh"
40 #include "G4NeutronHPDataUsed.hh"
41 
43  {
44  G4String tString = "/FS/";
45  G4bool dbool;
46  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
47  G4String filename = aFile.GetName();
48  SetAZMs( A, Z, M, aFile );
49  if(!dbool)
50  {
51  hasAnyData = false;
52  hasFSData = false;
53  hasXsec = false;
54  return;
55  }
56  std::ifstream theData(filename, std::ios::in);
57 
58  // here it comes
59  G4int infoType, dataType;
60  hasFSData = false;
61  while (theData >> infoType)
62  {
63  hasFSData = true;
64  theData >> dataType;
65  switch(infoType)
66  {
67  case 1:
68  if(dataType==4) theNeutronAngularDis.Init(theData);
69  if(dataType==5) thePromptNeutronEnDis.Init(theData);
70  if(dataType==12) theFinalStatePhotons.InitMean(theData);
71  if(dataType==14) theFinalStatePhotons.InitAngular(theData);
72  if(dataType==15) theFinalStatePhotons.InitEnergies(theData);
73  break;
74  case 2:
75  if(dataType==1) theFinalStateNeutrons.InitMean(theData);
76  break;
77  case 3:
78  if(dataType==1) theFinalStateNeutrons.InitDelayed(theData);
79  if(dataType==5) theDelayedNeutronEnDis.Init(theData);
80  break;
81  case 4:
82  if(dataType==1) theFinalStateNeutrons.InitPrompt(theData);
83  break;
84  case 5:
85  if(dataType==1) theEnergyRelease.Init(theData);
86  break;
87  default:
88  G4cout << "G4NeutronHPFSFissionFS::Init: unknown data type"<<dataType<<G4endl;
89  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPFSFissionFS::Init: unknown data type");
90  break;
91  }
92  }
93  targetMass = theFinalStateNeutrons.GetTargetMass();
94  theData.close();
95  }
96 
97 
99  G4int nDelayed, G4double * theDecayConst)
100  {
101  G4int i;
103  G4ReactionProduct boosted;
104  boosted.Lorentz(theNeutron, theTarget);
105  G4double eKinetic = boosted.GetKineticEnergy();
106 
107 // Build neutrons
108  G4ReactionProduct * theNeutrons = new G4ReactionProduct[nPrompt+nDelayed];
109  for(i=0; i<nPrompt+nDelayed; i++)
110  {
111  theNeutrons[i].SetDefinition(G4Neutron::Neutron());
112  }
113 
114 // sample energies
115  G4int it, dummy;
116  G4double tempE;
117  for(i=0; i<nPrompt; i++)
118  {
119  tempE = thePromptNeutronEnDis.Sample(eKinetic, dummy); // energy distribution (file5) always in lab
120  theNeutrons[i].SetKineticEnergy(tempE);
121  }
122  for(i=nPrompt; i<nPrompt+nDelayed; i++)
123  {
124  theNeutrons[i].SetKineticEnergy(theDelayedNeutronEnDis.Sample(eKinetic, it)); // dito
125  if(it==0) theNeutrons[i].SetKineticEnergy(thePromptNeutronEnDis.Sample(eKinetic, dummy));
126  theDecayConst[i-nPrompt] = theFinalStateNeutrons.GetDecayConstant(it); // this is returned
127  }
128 
129 // sample neutron angular distribution
130  for(i=0; i<nPrompt+nDelayed; i++)
131  {
132  theNeutronAngularDis.SampleAndUpdate(theNeutrons[i]); // angular comes back in lab automatically
133  }
134 
135 // already in lab. Add neutrons to dynamic particle vector
136  for(i=0; i<nPrompt+nDelayed; i++)
137  {
139  dp->SetDefinition(theNeutrons[i].GetDefinition());
140  dp->SetMomentum(theNeutrons[i].GetMomentum());
141  aResult->push_back(dp);
142  }
143  delete [] theNeutrons;
144 // return the result
145  return aResult;
146  }
147 
148 void G4NeutronHPFSFissionFS::SampleNeutronMult(G4int&all, G4int&Prompt, G4int&delayed, 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  {
157  Prompt = 0;
158  delayed = 0;
159  G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
160  all = G4Poisson(totalNeutronMulti-off);
161  all += off;
162  }
163  else
164  {
165  Prompt = G4Poisson(promptNeutronMulti-off);
166  Prompt += off;
167  delayed = G4Poisson(delayedNeutronMulti);
168  all = Prompt+delayed;
169  }
170 }
171 
173 {
174 // sample photons
176  G4ReactionProduct boosted;
177 // the photon distributions are in the Nucleus rest frame.
178  boosted.Lorentz(theNeutron, theTarget);
179  G4double anEnergy = boosted.GetKineticEnergy();
180  temp = theFinalStatePhotons.GetPhotons(anEnergy);
181  if(temp == 0) { return 0; }
182 
183 // lorentz transform, and add photons to final state
184  unsigned int i;
186  for(i=0; i<temp->size(); i++)
187  {
188  // back to lab
189  temp->operator[](i)->Lorentz(*(temp->operator[](i)), -1.*theTarget);
190  G4DynamicParticle * 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 }