Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPLabAngularEnergy.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// 080808 Bug fix in serching mu bin and index for theBuff2b by T. Koi
31//
34#include "G4SystemOfUnits.hh"
35#include "Randomize.hh"
36#include "G4Gamma.hh"
37#include "G4Electron.hh"
38#include "G4Positron.hh"
39#include "G4Neutron.hh"
40#include "G4Proton.hh"
41#include "G4Deuteron.hh"
42#include "G4Triton.hh"
43#include "G4He3.hh"
44#include "G4Alpha.hh"
45
46void G4NeutronHPLabAngularEnergy::Init(std::ifstream & aDataFile)
47{
48 aDataFile >> nEnergies;
49 theManager.Init(aDataFile);
50 theEnergies = new G4double[nEnergies];
51 nCosTh = new G4int[nEnergies];
52 theData = new G4NeutronHPVector * [nEnergies];
53 theSecondManager = new G4InterpolationManager [nEnergies];
54 for(G4int i=0; i<nEnergies; i++)
55 {
56 aDataFile >> theEnergies[i];
57 theEnergies[i]*=eV;
58 aDataFile >> nCosTh[i];
59 theSecondManager[i].Init(aDataFile);
60 theData[i] = new G4NeutronHPVector[nCosTh[i]];
61 G4double label;
62 for(G4int ii=0; ii<nCosTh[i]; ii++)
63 {
64 aDataFile >> label;
65 theData[i][ii].SetLabel(label);
66 theData[i][ii].Init(aDataFile, eV);
67 }
68 }
69}
70
72{
74 G4int Z = static_cast<G4int>(massCode/1000);
75 G4int A = static_cast<G4int>(massCode-1000*Z);
76
77 if(massCode==0)
78 {
80 }
81 else if(A==0)
82 {
84 if(Z==1) result->SetDefinition(G4Positron::Positron());
85 }
86 else if(A==1)
87 {
89 if(Z==1) result->SetDefinition(G4Proton::Proton());
90 }
91 else if(A==2)
92 {
94 }
95 else if(A==3)
96 {
98 if(Z==2) result->SetDefinition(G4He3::He3());
99 }
100 else if(A==4)
101 {
102 result->SetDefinition(G4Alpha::Alpha());
103 if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "Unknown ion case 1");
104 }
105 else
106 {
107 throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPLabAngularEnergy: Unknown ion case 2");
108 }
109
110 // get theta, E
111 G4double cosTh, secEnergy;
112 G4int i, it(0);
113 // find the energy bin
114 for(i=0; i<nEnergies; i++)
115 {
116 it = i;
117 if ( anEnergy < theEnergies[i] ) break;
118 }
119 //080808
120 //if ( it == 0 || it == nEnergies-1 ) // it marks the energy bin
121 if ( it == 0 ) // it marks the energy bin
122 {
123if(it==0) G4cout << "080808 Something unexpected is happen in G4NeutronHPLabAngularEnergy " << G4endl;
124 // integrate the prob for each costh, and select theta.
125 G4double * running = new G4double [nCosTh[it]];
126 running[0]=0;
127 for(i=0;i<nCosTh[it]; i++)
128 {
129 if(i!=0) running[i] = running[i-1];
130 running[i]+=theData[it][i].GetIntegral(); // Does interpolated integral.
131 }
132 G4double random = running[nCosTh[it]-1]*G4UniformRand();
133 G4int ith(0);
134 for(i=0;i<nCosTh[it]; i++)
135 {
136 ith = i;
137 if(random<running[i]) break;
138 }
139 //080807
140 //if ( ith == 0 || ith == nCosTh[it]-1 ) //ith marks the angluar bin
141 if ( ith == 0 ) //ith marks the angluar bin
142 {
143 cosTh = theData[it][ith].GetLabel();
144 secEnergy = theData[it][ith].Sample();
145 currentMeanEnergy = theData[it][ith].GetMeanX();
146 }
147 else
148 {
149 //080808
150 //G4double x1 = theData[it][ith-1].GetIntegral();
151 //G4double x2 = theData[it][ith].GetIntegral();
152 G4double x1 = running [ ith-1 ];
153 G4double x2 = running [ ith ];
154 G4double x = random;
155 G4double y1 = theData[it][ith-1].GetLabel();
156 G4double y2 = theData[it][ith].GetLabel();
157 cosTh = theInt.Interpolate(theSecondManager[it].GetInverseScheme(ith),
158 x, x1, x2, y1, y2);
159 G4NeutronHPVector theBuff1;
160 theBuff1.SetInterpolationManager(theData[it][ith-1].GetInterpolationManager());
161 G4NeutronHPVector theBuff2;
162 theBuff2.SetInterpolationManager(theData[it][ith].GetInterpolationManager());
163 x1=y1;
164 x2=y2;
165 G4double y, mu;
166 for(i=0;i<theData[it][ith-1].GetVectorLength(); i++)
167 {
168 mu = theData[it][ith-1].GetX(i);
169 y1 = theData[it][ith-1].GetY(i);
170 y2 = theData[it][ith].GetY(mu);
171
172 y = theInt.Interpolate(theSecondManager[it].GetScheme(ith),
173 cosTh, x1,x2,y1,y2);
174 theBuff1.SetData(i, mu, y);
175 }
176 for(i=0;i<theData[it][ith].GetVectorLength(); i++)
177 {
178 mu = theData[it][ith].GetX(i);
179 y1 = theData[it][ith-1].GetY(mu);
180 y2 = theData[it][ith].GetY(i);
181 y = theInt.Interpolate(theSecondManager[it].GetScheme(ith),
182 cosTh, x1,x2,y1,y2);
183 theBuff2.SetData(i, mu, y);
184 }
185 G4NeutronHPVector theStore;
186 theStore.Merge(&theBuff1, &theBuff2);
187 secEnergy = theStore.Sample();
188 currentMeanEnergy = theStore.GetMeanX();
189 }
190 delete [] running;
191 }
192 else // this is the small big else.
193 {
194 G4double x, x1, x2, y1, y2, y, tmp, E;
195 // integrate the prob for each costh, and select theta.
197 run1.SetY(0, 0.);
198 for(i=0;i<nCosTh[it-1]; i++)
199 {
200 if(i!=0) run1.SetY(i, run1.GetY(i-1));
201 run1.SetX(i, theData[it-1][i].GetLabel());
202 run1.SetY(i, run1.GetY(i)+theData[it-1][i].GetIntegral());
203 }
205 run2.SetY(0, 0.);
206 for(i=0;i<nCosTh[it]; i++)
207 {
208 if(i!=0) run2.SetY(i, run2.GetY(i-1));
209 run2.SetX(i, theData[it][i].GetLabel());
210 run2.SetY(i, run2.GetY(i)+theData[it][i].GetIntegral());
211 }
212 // get the distributions for the correct neutron energy
213 x = anEnergy;
214 x1 = theEnergies[it-1];
215 x2 = theEnergies[it];
216 G4NeutronHPVector thBuff1; // to be interpolated as run1.
217 thBuff1.SetInterpolationManager(theSecondManager[it-1]);
218 for(i=0; i<run1.GetVectorLength(); i++)
219 {
220 tmp = run1.GetX(i); //theta
221 y1 = run1.GetY(i); // integral
222 y2 = run2.GetY(tmp);
223 y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
224 thBuff1.SetData(i, tmp, y);
225 }
226 G4NeutronHPVector thBuff2;
227 thBuff2.SetInterpolationManager(theSecondManager[it]);
228 for(i=0; i<run2.GetVectorLength(); i++)
229 {
230 tmp = run2.GetX(i); //theta
231 y1 = run1.GetY(tmp); // integral
232 y2 = run2.GetY(i);
233 y = theInt.Lin(x, x1,x2,y1,y2);
234 thBuff2.SetData(i, tmp, y);
235 }
236 G4NeutronHPVector theThVec;
237 theThVec.Merge(&thBuff1 ,&thBuff2); // takes care of interpolation
238 G4double random = (theThVec.GetY(theThVec.GetVectorLength()-1)
239 -theThVec.GetY(0)) *G4UniformRand();
240 G4int ith(0);
241 for(i=1;i<theThVec.GetVectorLength(); i++)
242 {
243 ith = i;
244 if(random<theThVec.GetY(i)-theThVec.GetY(0)) break;
245 }
246 {
247 // calculate theta
248 G4double xx, xx1, xx2, yy1, yy2;
249 xx = random;
250 xx1 = theThVec.GetY(ith-1)-theThVec.GetY(0); // integrals
251 xx2 = theThVec.GetY(ith)-theThVec.GetY(0);
252 yy1 = theThVec.GetX(ith-1); // std::cos(theta)
253 yy2 = theThVec.GetX(ith);
254 cosTh = theInt.Interpolate(theSecondManager[it].GetScheme(ith),
255 xx, xx1,xx2,yy1,yy2);
256 }
257 G4int i1(0), i2(0);
258 // get the indixes of the vectors close to theta for low energy
259 // first it-1 !!!! i.e. low in energy
260 for(i=0; i<nCosTh[it-1]; i++)
261 {
262 i1 = i;
263 if(cosTh<theData[it-1][i].GetLabel()) break;
264 }
265 // now get the prob at this energy for the right theta value
266 x = cosTh;
267 x1 = theData[it-1][i1-1].GetLabel();
268 x2 = theData[it-1][i1].GetLabel();
269 G4NeutronHPVector theBuff1a;
270 theBuff1a.SetInterpolationManager(theData[it-1][i1-1].GetInterpolationManager());
271 for(i=0;i<theData[it-1][i1-1].GetVectorLength(); i++)
272 {
273 E = theData[it-1][i1-1].GetX(i);
274 y1 = theData[it-1][i1-1].GetY(i);
275 y2 = theData[it-1][i1].GetY(E);
276 y = theInt.Lin(x, x1,x2,y1,y2);
277 theBuff1a.SetData(i, E, y); // wrong E, right theta.
278 }
279 G4NeutronHPVector theBuff2a;
280 theBuff2a.SetInterpolationManager(theData[it-1][i1].GetInterpolationManager());
281 for(i=0;i<theData[it-1][i1].GetVectorLength(); i++)
282 {
283 E = theData[it-1][i1].GetX(i);
284 y1 = theData[it-1][i1-1].GetY(E);
285 y2 = theData[it-1][i1].GetY(i);
286 y = theInt.Lin(x, x1,x2,y1,y2);
287 theBuff2a.SetData(i, E, y); // wrong E, right theta.
288 }
289 G4NeutronHPVector theStore1;
290 theStore1.Merge(&theBuff1a, &theBuff2a); // wrong E, right theta, complete binning
291
292 // get the indixes of the vectors close to theta for high energy
293 // then it !!!! i.e. high in energy
294 for(i=0; i<nCosTh[it]; i++)
295 {
296 i2 = i;
297 if(cosTh<theData[it][i2].GetLabel()) break;
298 } // sonderfaelle mit i1 oder i2 head on fehlen. @@@@@
299 x1 = theData[it][i2-1].GetLabel();
300 x2 = theData[it][i2].GetLabel();
301 G4NeutronHPVector theBuff1b;
302 theBuff1b.SetInterpolationManager(theData[it][i2-1].GetInterpolationManager());
303 for(i=0;i<theData[it][i2-1].GetVectorLength(); i++)
304 {
305 E = theData[it][i2-1].GetX(i);
306 y1 = theData[it][i2-1].GetY(i);
307 y2 = theData[it][i2].GetY(E);
308 y = theInt.Lin(x, x1,x2,y1,y2);
309 theBuff1b.SetData(i, E, y); // wrong E, right theta.
310 }
311 G4NeutronHPVector theBuff2b;
312 theBuff2b.SetInterpolationManager(theData[it][i2].GetInterpolationManager());
313 //080808 i1 -> i2
314 //for(i=0;i<theData[it][i1].GetVectorLength(); i++)
315 for(i=0;i<theData[it][i2].GetVectorLength(); i++)
316 {
317 //E = theData[it][i1].GetX(i);
318 //y1 = theData[it][i1-1].GetY(E);
319 //y2 = theData[it][i1].GetY(i);
320 E = theData[it][i2].GetX(i);
321 y1 = theData[it][i2-1].GetY(E);
322 y2 = theData[it][i2].GetY(i);
323 y = theInt.Lin(x, x1,x2,y1,y2);
324 theBuff2b.SetData(i, E, y); // wrong E, right theta.
325 }
326 G4NeutronHPVector theStore2;
327 theStore2.Merge(&theBuff1b, &theBuff2b); // wrong E, right theta, complete binning
328 // now get to the right energy.
329
330 x = anEnergy;
331 x1 = theEnergies[it-1];
332 x2 = theEnergies[it];
333 G4NeutronHPVector theOne1;
335 for(i=0; i<theStore1.GetVectorLength(); i++)
336 {
337 E = theStore1.GetX(i);
338 y1 = theStore1.GetY(i);
339 y2 = theStore2.GetY(E);
340 y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
341 theOne1.SetData(i, E, y); // both correct
342 }
343 G4NeutronHPVector theOne2;
345 for(i=0; i<theStore2.GetVectorLength(); i++)
346 {
347 E = theStore2.GetX(i);
348 y1 = theStore1.GetY(E);
349 y2 = theStore2.GetY(i);
350 y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
351 theOne2.SetData(i, E, y); // both correct
352 }
353 G4NeutronHPVector theOne;
354 theOne.Merge(&theOne1, &theOne2); // both correct, complete binning
355
356 secEnergy = theOne.Sample();
357 currentMeanEnergy = theOne.GetMeanX();
358 }
359
360// now do random direction in phi, and fill the result.
361
362 result->SetKineticEnergy(secEnergy);
363
364 G4double phi = twopi*G4UniformRand();
365 G4double theta = std::acos(cosTh);
366 G4double sinth = std::sin(theta);
367 G4double mtot = result->GetTotalMomentum();
368 G4ThreeVector tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) );
369 result->SetMomentum(tempVector);
370
371 return result;
372}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4He3 * He3()
Definition: G4He3.cc:94
void Init(G4int aScheme, G4int aRange)
G4InterpolationScheme GetScheme(G4int index) const
G4double Lin(G4double x, G4double x1, G4double x2, G4double y1, G4double y2)
G4double Interpolate(G4InterpolationScheme aScheme, G4double x, G4double x1, G4double x2, G4double y1, G4double y2) const
void Init(std::ifstream &aDataFile)
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass)
void Init(std::ifstream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
void SetLabel(G4double aLabel)
void Merge(G4NeutronHPVector *active, G4NeutronHPVector *passive)
void SetX(G4int i, G4double e)
G4int GetVectorLength() const
G4double GetX(G4int i) const
const G4InterpolationManager & GetInterpolationManager() const
G4double GetY(G4double x)
void SetData(G4int i, G4double x, G4double y)
void SetInterpolationManager(const G4InterpolationManager &aManager)
void SetY(G4int i, G4double x)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetTotalMomentum() const
void SetKineticEnergy(const G4double en)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
static G4Triton * Triton()
Definition: G4Triton.cc:95