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