Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPContAngularPar.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 // 09-May-06 fix in Sample by T. Koi
31 // 080318 Fix Compilation warnings - gcc-4.3.0 by T. Koi
32 // (This fix has a real effect to the code.)
33 // 080409 Fix div0 error with G4FPE by T. Koi
34 // 080612 Fix contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #1
35 // 080714 Limiting the sum of energy of secondary particles by T. Koi
36 // 080801 Fix div0 error wiht G4FPE and memory leak by T. Koi
37 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
38 //
39 
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
44 #include "G4Gamma.hh"
45 #include "G4Electron.hh"
46 #include "G4Positron.hh"
47 #include "G4Neutron.hh"
48 #include "G4Proton.hh"
49 #include "G4Deuteron.hh"
50 #include "G4Triton.hh"
51 #include "G4He3.hh"
52 #include "G4Alpha.hh"
53 #include "G4NeutronHPVector.hh"
54 #include "G4NucleiProperties.hh"
56 #include "G4IonTable.hh"
57 
58  void G4NeutronHPContAngularPar::Init(std::istream & aDataFile)
59  {
60  aDataFile >> theEnergy >> nEnergies >> nDiscreteEnergies >> nAngularParameters;
61  theEnergy *= eV;
62  theAngular = new G4NeutronHPList [nEnergies];
63  for(G4int i=0; i<nEnergies; i++)
64  {
65  G4double sEnergy;
66  aDataFile >> sEnergy;
67  sEnergy*=eV;
68  theAngular[i].SetLabel(sEnergy);
69  theAngular[i].Init(aDataFile, nAngularParameters, 1.);
70  }
71  }
72 
74  G4NeutronHPContAngularPar::Sample(G4double anEnergy, G4double massCode, G4double /*targetMass*/,
75  G4int angularRep, G4int /*interpolE*/ )
76  {
77  G4ReactionProduct * result = new G4ReactionProduct;
78  G4int Z = static_cast<G4int>(massCode/1000);
79  G4int A = static_cast<G4int>(massCode-1000*Z);
80  if(massCode==0)
81  {
82  result->SetDefinition(G4Gamma::Gamma());
83  }
84  else if(A==0)
85  {
87  if(Z==1) result->SetDefinition(G4Positron::Positron());
88  }
89  else if(A==1)
90  {
92  if(Z==1) result->SetDefinition(G4Proton::Proton());
93  }
94  else if(A==2)
95  {
97  }
98  else if(A==3)
99  {
100  result->SetDefinition(G4Triton::Triton());
101  if(Z==2) result->SetDefinition(G4He3::He3());
102  }
103  else if(A==4)
104  {
105  result->SetDefinition(G4Alpha::Alpha());
106  if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPContAngularPar: Unknown ion case 1");
107  }
108  else
109  {
110  //result->SetDefinition(G4ParticleTable::GetParticleTable()->FindIon(Z,A,0,Z));
111  result->SetDefinition(G4IonTable::GetIonTable()->GetIon(Z,A,0));
112  }
113  G4int i(0);
114  G4int it(0);
115  G4double fsEnergy(0);
116  G4double cosTh(0);
117 
118  if( angularRep == 1 )
119  {
120 // 080612 Fix contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #1
121  //if (interpolE == 2)
122 //110609 above was wrong interupition, pointed out by E.Mendoza and D.Cano (CIMAT)
123 //Following are reviesd version written by T.Koi (SLAC)
124  if ( nDiscreteEnergies != 0 )
125  {
126 
127 //1st check remaining_energy
128 // if this is the first set it. (How?)
129  if ( fresh == true )
130  {
131  //Discrete Lines, larger energies come first
132  //Continues Emssions, low to high LAST
133  remaining_energy = std::max ( theAngular[0].GetLabel() , theAngular[nEnergies-1].GetLabel() );
134  fresh = false;
135  }
136 
137  //Cheating for small remaining_energy
138  //TEMPORAL SOLUTION
139  if ( nDiscreteEnergies == nEnergies )
140  {
141  remaining_energy = std::max ( remaining_energy , theAngular[nDiscreteEnergies-1].GetLabel() ); //Minimum Line
142  }
143  else
144  {
145  //G4double cont_min = theAngular[nDiscreteEnergies].GetLabel();
146  //if ( theAngular[nDiscreteEnergies].GetLabel() == 0.0 ) cont_min = theAngular[nDiscreteEnergies+1].GetLabel();
147  G4double cont_min=0.0;
148  for ( G4int j = nDiscreteEnergies ; j < nEnergies ; j++ )
149  {
150  cont_min = theAngular[j].GetLabel();
151  if ( theAngular[j].GetValue(0) != 0.0 ) break;
152  }
153  remaining_energy = std::max ( remaining_energy , std::min ( theAngular[nDiscreteEnergies-1].GetLabel() , cont_min ) ); //Minimum Line or grid
154  }
155 //
156  G4double random = G4UniformRand();
157 
158  G4double * running = new G4double[nEnergies+1];
159  running[0] = 0.0;
160 
161  for ( G4int j = 0 ; j < nDiscreteEnergies ; j++ )
162  {
163  G4double delta = 0.0;
164  if ( theAngular[j].GetLabel() <= remaining_energy ) delta = theAngular[i].GetValue(0);
165  running[j+1] = running[j] + delta;
166  }
167  G4double tot_prob_DIS = running[ nDiscreteEnergies ];
168 
169  for ( G4int j = nDiscreteEnergies ; j < nEnergies ; j++ )
170  {
171  G4double delta = 0.0;
172  G4double e_low = 0.0;
173  G4double e_high = 0.0;
174  if ( theAngular[j].GetLabel() <= remaining_energy ) delta = theAngular[j].GetValue(0);
175 
176  //To calculate Prob. e_low and e_high should be in eV
177  //There are two case
178  //1:theAngular[nDiscreteEnergies].GetLabel() != 0.0
179  // delta should be used between j-1 and j
180  // At j = nDiscreteEnergies (the first) e_low should be set explicitly
181  if ( theAngular[j].GetLabel() != 0 )
182  {
183  if ( j == nDiscreteEnergies ) {
184  e_low = 0.0/eV;
185  } else {
186  e_low = theAngular[j-1].GetLabel()/eV;
187  }
188  e_high = theAngular[j].GetLabel()/eV;
189  }
190  //2:theAngular[nDiscreteEnergies].GetLabel() == 0.0
191  // delta should be used between j and j+1
192  if ( theAngular[j].GetLabel() == 0.0 ) {
193  e_low = theAngular[j].GetLabel()/eV;
194  if ( j != nEnergies-1 ) {
195  e_high = theAngular[j+1].GetLabel()/eV;
196  } else {
197  e_high = theAngular[j].GetLabel()/eV;
198  if ( theAngular[j].GetValue(0) != 0.0 ) {
199  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPContAngularPar: Unexpected non zero value of theAngular[nEnergies-1].GetValue(0)");
200  }
201  }
202  }
203 
204  running[j+1] = running[j] + ( ( e_high - e_low ) * delta );
205  }
206  G4double tot_prob_CON = running[ nEnergies ] - running[ nDiscreteEnergies ];
207 
208 /*
209  For FPE debugging
210  if (tot_prob_DIS + tot_prob_CON == 0 ) {
211  G4cout << "TKDB tot_prob_DIS + tot_prob_CON " << tot_prob_DIS + tot_prob_CON << G4endl;
212  G4cout << "massCode " << massCode << G4endl;
213  G4cout << "nDiscreteEnergies " << nDiscreteEnergies << " nEnergies " << nEnergies << G4endl;
214  for ( int j = nDiscreteEnergies ; j < nEnergies ; j++ ) {
215  G4cout << j << " " << theAngular[j].GetLabel() << " " << theAngular[j].GetValue(0) << G4endl;
216  }
217  }
218 */
219  // Normalize random
220  random *= (tot_prob_DIS + tot_prob_CON);
221 //2nd Judge Discrete or not This shoudl be relatively close to 1 For safty
222  if ( random <= ( tot_prob_DIS / ( tot_prob_DIS + tot_prob_CON ) ) || nDiscreteEnergies == nEnergies )
223  {
224 // Discrete Emission
225  for ( G4int j = 0 ; j < nDiscreteEnergies ; j++ )
226  {
227  //Here we should use i+1
228  if ( random < running[ j+1 ] )
229  {
230  it = j;
231  break;
232  }
233  }
234  fsEnergy = theAngular[ it ].GetLabel();
235 
236  G4NeutronHPLegendreStore theStore(1);
237  theStore.Init(0,fsEnergy,nAngularParameters);
238  for (G4int j=0;j<nAngularParameters;j++)
239  {
240  theStore.SetCoeff(0,j,theAngular[it].GetValue(j));
241  }
242  // use it to sample.
243  cosTh = theStore.SampleMax(fsEnergy);
244  //Done
245  }
246  else
247  {
248 // Continuous Emission
249  for ( G4int j = nDiscreteEnergies ; j < nEnergies ; j++ )
250  {
251  //Here we should use i
252  if ( random < running[ j ] )
253  {
254  it = j;
255  break;
256  }
257  }
258 
259  G4double x1 = running[it-1];
260  G4double x2 = running[it];
261 
262  G4double y1 = 0.0;
263  if ( it != nDiscreteEnergies )
264  y1 = theAngular[it-1].GetLabel();
265  G4double y2 = theAngular[it].GetLabel();
266 
267  fsEnergy = theInt.Interpolate(theManager.GetInverseScheme(it),
268  random,x1,x2,y1,y2);
269 
270  G4NeutronHPLegendreStore theStore(2);
271  theStore.Init(0,y1,nAngularParameters);
272  theStore.Init(1,y2,nAngularParameters);
273  theStore.SetManager(theManager);
274  for (G4int j=0;j<nAngularParameters;j++)
275  {
276  G4int itt = it;
277  if ( it == nDiscreteEnergies ) itt = it+1; //"This case "it-1" has data for Discrete, so we will use an extrpolate values it and it+1
278  if ( it == 0 )
279  {
280  //Safty for unexpected it = 0;
281  //G4cout << "110611 G4NeutronHPContAngularPar::Sample it = 0; invetigation required " << G4endl;
282  itt = it+1;
283  }
284  theStore.SetCoeff(0,j,theAngular[itt-1].GetValue(j));
285  theStore.SetCoeff(1,j,theAngular[itt].GetValue(j));
286  }
287  // use it to sample.
288  cosTh = theStore.SampleMax(fsEnergy);
289 
290  //Done
291  }
292 
293  //TK080711
294  remaining_energy -= fsEnergy;
295  //TK080711
296 
297  //080801b
298  delete[] running;
299  //080801b
300  }
301  else
302  {
303  // Only continue, TK will clean up
304 
305  //080714
306  if ( fresh == true )
307  {
308  remaining_energy = theAngular[ nEnergies-1 ].GetLabel();
309  fresh = false;
310  }
311  //080714
312  G4double random = G4UniformRand();
313  G4double * running = new G4double[nEnergies];
314  running[0]=0;
315  G4double weighted = 0;
316  for(i=1; i<nEnergies; i++)
317  {
318 /*
319  if(i!=0)
320  {
321  running[i]=running[i-1];
322  }
323  running[i] += theInt.GetBinIntegral(theManager.GetScheme(i-1),
324  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
325  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
326  weighted += theInt.GetWeightedBinIntegral(theManager.GetScheme(i-1),
327  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
328  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
329 */
330 
331  running[i]=running[i-1];
332  if ( remaining_energy >= theAngular[i].GetLabel() )
333  {
334  running[i] += theInt.GetBinIntegral(theManager.GetScheme(i-1),
335  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
336  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
337  weighted += theInt.GetWeightedBinIntegral(theManager.GetScheme(i-1),
338  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
339  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
340  }
341  }
342  // cash the mean energy in this distribution
343  //080409 TKDB
344  if ( nEnergies == 1 || running[nEnergies-1] == 0 )
345  currentMeanEnergy = 0.0;
346  else
347  {
348  currentMeanEnergy = weighted/running[nEnergies-1];
349  }
350 
351  //080409 TKDB
352  if ( nEnergies == 1 ) it = 0;
353 
354  //080729
355  if ( running[nEnergies-1] != 0 )
356  {
357  for ( i = 1 ; i < nEnergies ; i++ )
358  {
359  it = i;
360  if ( random < running [ i ] / running [ nEnergies-1 ] ) break;
361  }
362  }
363 
364  //080714
365  if ( running [ nEnergies-1 ] == 0 ) it = 0;
366  //080714
367 
368  if (it<nDiscreteEnergies||it==0)
369  {
370  if(it == 0)
371  {
372  fsEnergy = theAngular[it].GetLabel();
373  G4NeutronHPLegendreStore theStore(1);
374  theStore.Init(0,fsEnergy,nAngularParameters);
375  for(i=0;i<nAngularParameters;i++)
376  {
377  theStore.SetCoeff(0,i,theAngular[it].GetValue(i));
378  }
379  // use it to sample.
380  cosTh = theStore.SampleMax(fsEnergy);
381  }
382  else
383  {
384  G4double e1, e2;
385  e1 = theAngular[it-1].GetLabel();
386  e2 = theAngular[it].GetLabel();
387  fsEnergy = theInt.Interpolate(theManager.GetInverseScheme(it),
388  random,
389  running[it-1]/running[nEnergies-1],
390  running[it]/running[nEnergies-1],
391  e1, e2);
392  // fill a Legendrestore
393  G4NeutronHPLegendreStore theStore(2);
394  theStore.Init(0,e1,nAngularParameters);
395  theStore.Init(1,e2,nAngularParameters);
396  for(i=0;i<nAngularParameters;i++)
397  {
398  theStore.SetCoeff(0,i,theAngular[it-1].GetValue(i));
399  theStore.SetCoeff(1,i,theAngular[it].GetValue(i));
400  }
401  // use it to sample.
402  theStore.SetManager(theManager);
403  cosTh = theStore.SampleMax(fsEnergy);
404  }
405  }
406  else // continuum contribution
407  {
408  G4double x1 = running[it-1]/running[nEnergies-1];
409  G4double x2 = running[it]/running[nEnergies-1];
410  G4double y1 = theAngular[it-1].GetLabel();
411  G4double y2 = theAngular[it].GetLabel();
412  fsEnergy = theInt.Interpolate(theManager.GetInverseScheme(it),
413  random,x1,x2,y1,y2);
414  G4NeutronHPLegendreStore theStore(2);
415  theStore.Init(0,y1,nAngularParameters);
416  theStore.Init(1,y2,nAngularParameters);
417  theStore.SetManager(theManager);
418  for(i=0;i<nAngularParameters;i++)
419  {
420  theStore.SetCoeff(0,i,theAngular[it-1].GetValue(i));
421  theStore.SetCoeff(1,i,theAngular[it].GetValue(i));
422  }
423  // use it to sample.
424  cosTh = theStore.SampleMax(fsEnergy);
425  }
426  delete [] running;
427 
428  //080714
429  remaining_energy -= fsEnergy;
430  //080714
431  }
432  }
433  else if(angularRep==2)
434  {
435  // first get the energy (already the right for this incoming energy)
436  G4int j;
437  G4double * running = new G4double[nEnergies];
438  running[0]=0;
439  G4double weighted = 0;
440  for(j=1; j<nEnergies; j++)
441  {
442  if(j!=0) running[j]=running[j-1];
443  running[j] += theInt.GetBinIntegral(theManager.GetScheme(j-1),
444  theAngular[j-1].GetLabel(), theAngular[j].GetLabel(),
445  theAngular[j-1].GetValue(0), theAngular[j].GetValue(0));
446  weighted += theInt.GetWeightedBinIntegral(theManager.GetScheme(j-1),
447  theAngular[j-1].GetLabel(), theAngular[j].GetLabel(),
448  theAngular[j-1].GetValue(0), theAngular[j].GetValue(0));
449  }
450  // cash the mean energy in this distribution
451  //080409 TKDB
452  //currentMeanEnergy = weighted/running[nEnergies-1];
453  if ( nEnergies == 1 )
454  currentMeanEnergy = 0.0;
455  else
456  currentMeanEnergy = weighted/running[nEnergies-1];
457 
458  G4int itt(0);
459  G4double randkal = G4UniformRand();
460  //080409 TKDB
461  //for(i=0; i<nEnergies; i++)
462  for(j=1; j<nEnergies; j++)
463  {
464  itt = j;
465  if(randkal<running[j]/running[nEnergies-1]) break;
466  }
467 
468  // interpolate the secondary energy.
469  G4double x, x1,x2,y1,y2;
470  if(itt==0) itt=1;
471  x = randkal*running[nEnergies-1];
472  x1 = running[itt-1];
473  x2 = running[itt];
474  G4double compoundFraction;
475  // interpolate energy
476  y1 = theAngular[itt-1].GetLabel();
477  y2 = theAngular[itt].GetLabel();
478  fsEnergy = theInt.Interpolate(theManager.GetInverseScheme(itt-1),
479  x, x1,x2,y1,y2);
480  // for theta interpolate the compoundFractions
481  G4double cLow = theAngular[itt-1].GetValue(1);
482  G4double cHigh = theAngular[itt].GetValue(1);
483  compoundFraction = theInt.Interpolate(theManager.GetScheme(itt),
484  fsEnergy, y1, y2, cLow,cHigh);
485  delete [] running;
486 
487  // get cosTh
488  G4double incidentEnergy = anEnergy;
489  G4double incidentMass = G4Neutron::Neutron()->GetPDGMass();
490  G4double productEnergy = fsEnergy;
491  G4double productMass = result->GetMass();
492  G4int targetZ = G4int(theTargetCode/1000);
493  G4int targetA = G4int(theTargetCode-1000*targetZ);
494  // To correspond to natural composition (-nat-) data files.
495  if ( targetA == 0 )
496  targetA = G4int ( theTarget->GetMass()/amu_c2 + 0.5 );
497  G4double targetMass = theTarget->GetMass();
498  G4int residualA = targetA+1-A;
499  G4int residualZ = targetZ-Z;
500  G4double residualMass = residualZ*G4Proton::Proton()->GetPDGMass();
501  residualMass +=(residualA-residualZ)*G4Neutron::Neutron()->GetPDGMass();
502  residualMass -= G4NucleiProperties::GetBindingEnergy( residualA , residualZ );
503  G4NeutronHPKallbachMannSyst theKallbach(compoundFraction,
504  incidentEnergy, incidentMass,
505  productEnergy, productMass,
506  residualMass, residualA, residualZ,
507  targetMass, targetA, targetZ);
508  cosTh = theKallbach.Sample(anEnergy);
509  }
510  else if(angularRep>10&&angularRep<16)
511  {
512  G4double random = G4UniformRand();
513  G4double * running = new G4double[nEnergies];
514  running[0]=0;
515  G4double weighted = 0;
516  for(i=1; i<nEnergies; i++)
517  {
518  if(i!=0) running[i]=running[i-1];
519  running[i] += theInt.GetBinIntegral(theManager.GetScheme(i-1),
520  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
521  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
522  weighted += theInt.GetWeightedBinIntegral(theManager.GetScheme(i-1),
523  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
524  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
525  }
526  // cash the mean energy in this distribution
527  //currentMeanEnergy = weighted/running[nEnergies-1];
528  if ( nEnergies == 1 )
529  currentMeanEnergy = 0.0;
530  else
531  currentMeanEnergy = weighted/running[nEnergies-1];
532 
533  //080409 TKDB
534  if ( nEnergies == 1 ) it = 0;
535  //for(i=0; i<nEnergies; i++)
536  for(i=1; i<nEnergies; i++)
537  {
538  it = i;
539  if(random<running[i]/running[nEnergies-1]) break;
540  }
541  if(it<nDiscreteEnergies||it==0)
542  {
543  if(it==0)
544  {
545  fsEnergy = theAngular[0].GetLabel();
546  G4NeutronHPVector theStore;
547  G4int aCounter = 0;
548  for(G4int j=1; j<nAngularParameters; j+=2)
549  {
550  theStore.SetX(aCounter, theAngular[0].GetValue(j));
551  theStore.SetY(aCounter, theAngular[0].GetValue(j+1));
552  aCounter++;
553  }
555  aMan.Init(angularRep-10, nAngularParameters-1);
556  theStore.SetInterpolationManager(aMan);
557  cosTh = theStore.Sample();
558  }
559  else
560  {
561  fsEnergy = theAngular[it].GetLabel();
562  G4NeutronHPVector theStore;
564  aMan.Init(angularRep-10, nAngularParameters-1);
565  theStore.SetInterpolationManager(aMan); // Store interpolates f(costh)
566  G4InterpolationScheme currentScheme = theManager.GetInverseScheme(it);
567  G4int aCounter = 0;
568  for(G4int j=1; j<nAngularParameters; j+=2)
569  {
570  theStore.SetX(aCounter, theAngular[it].GetValue(j));
571  theStore.SetY(aCounter, theInt.Interpolate(currentScheme,
572  random,
573  running[it-1]/running[nEnergies-1],
574  running[it]/running[nEnergies-1],
575  theAngular[it-1].GetValue(j+1),
576  theAngular[it].GetValue(j+1)));
577  aCounter++;
578  }
579  cosTh = theStore.Sample();
580  }
581  }
582  else
583  {
584  G4double x1 = running[it-1]/running[nEnergies-1];
585  G4double x2 = running[it]/running[nEnergies-1];
586  G4double y1 = theAngular[it-1].GetLabel();
587  G4double y2 = theAngular[it].GetLabel();
588  fsEnergy = theInt.Interpolate(theManager.GetInverseScheme(it),
589  random,x1,x2,y1,y2);
590  G4NeutronHPVector theBuff1;
591  G4NeutronHPVector theBuff2;
593  aMan.Init(angularRep-10, nAngularParameters-1);
594 // theBuff1.SetInterpolationManager(aMan); // Store interpolates f(costh)
595 // theBuff2.SetInterpolationManager(aMan); // Store interpolates f(costh)
596 // Bug Report #1366 from L. Russell
597  //for(i=0; i<nAngularParameters; i++) // i=1 ist wichtig!
598  //{
599  // theBuff1.SetX(i, theAngular[it-1].GetValue(i));
600  // theBuff1.SetY(i, theAngular[it-1].GetValue(i+1));
601  // theBuff2.SetX(i, theAngular[it].GetValue(i));
602  // theBuff2.SetY(i, theAngular[it].GetValue(i+1));
603  // i++;
604  //}
605  {
606  G4int j;
607  for(i=0,j=1; i<nAngularParameters; i++,j+=2)
608  {
609  theBuff1.SetX(i, theAngular[it-1].GetValue(j));
610  theBuff1.SetY(i, theAngular[it-1].GetValue(j+1));
611  theBuff2.SetX(i, theAngular[it].GetValue(j));
612  theBuff2.SetY(i, theAngular[it].GetValue(j+1));
613  }
614  }
615  G4NeutronHPVector theStore;
616  theStore.SetInterpolationManager(aMan); // Store interpolates f(costh)
617  x1 = y1;
618  x2 = y2;
619  G4double x, y;
620  //for(i=0;i<theBuff1.GetVectorLength(); i++);
621  for(i=0;i<theBuff1.GetVectorLength(); i++)
622  {
623  x = theBuff1.GetX(i); // costh binning identical
624  y1 = theBuff1.GetY(i);
625  y2 = theBuff2.GetY(i);
626  y = theInt.Interpolate(theManager.GetScheme(it),
627  fsEnergy, theAngular[it-1].GetLabel(),
628  theAngular[it].GetLabel(), y1, y2);
629  theStore.SetX(i, x);
630  theStore.SetY(i, y);
631  }
632  cosTh = theStore.Sample();
633  }
634  delete [] running;
635  }
636  else
637  {
638  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPContAngularPar::Sample: Unknown angular representation");
639  }
640  result->SetKineticEnergy(fsEnergy);
641  G4double phi = twopi*G4UniformRand();
642  G4double theta = std::acos(cosTh);
643  G4double sinth = std::sin(theta);
644  G4double mtot = result->GetTotalMomentum();
645  G4ThreeVector tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) );
646  result->SetMomentum(tempVector);
647 // return the result.
648  return result;
649  }
G4double GetY(G4double x)
G4int GetVectorLength() const
void SetInterpolationManager(const G4InterpolationManager &aManager)
G4double SampleMax(G4double energy)
G4double GetTotalMomentum() const
G4double GetLabel()
void SetLabel(G4double aLabel)
void Init(G4int aScheme, G4int aRange)
void SetKineticEnergy(const G4double en)
void Init(G4int i, G4double e, G4int n)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetCoeff(G4int i, G4int l, G4double coeff)
G4double GetX(G4int i) const
int G4int
Definition: G4Types.hh:78
G4double Interpolate(G4InterpolationScheme aScheme, G4double x, G4double x1, G4double x2, G4double y1, G4double y2) const
void SetY(G4int i, G4double x)
G4double GetWeightedBinIntegral(const G4InterpolationScheme &aScheme, const G4double x1, const G4double x2, const G4double y1, const G4double y2)
#define G4UniformRand()
Definition: Randomize.hh:87
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass, G4int angularRep, G4int interpol)
G4InterpolationScheme GetScheme(G4int index) const
G4double GetBinIntegral(const G4InterpolationScheme &aScheme, const G4double x1, const G4double x2, const G4double y1, const G4double y2)
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetX(G4int i, G4double e)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:80
G4double Sample(G4double anEnergy)
void Init(std::istream &aDataFile, G4int nPar, G4double unit=1.)
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4InterpolationScheme
G4double GetPDGMass() const
static G4double GetBindingEnergy(const G4int A, const G4int Z)
void Init(std::istream &aDataFile)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetManager(G4InterpolationManager &aManager)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GetValue(G4int i)
G4InterpolationScheme GetInverseScheme(G4int index)
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
double G4double
Definition: G4Types.hh:76
float amu_c2
Definition: hepunit.py:277
static G4He3 * He3()
Definition: G4He3.cc:94
G4double GetMass() const