Geant4  10.02
G4ParticleHPContAngularPar.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 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
40 //
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
45 #include "G4Gamma.hh"
46 #include "G4Electron.hh"
47 #include "G4Positron.hh"
48 #include "G4Neutron.hh"
49 #include "G4Proton.hh"
50 #include "G4Deuteron.hh"
51 #include "G4Triton.hh"
52 #include "G4He3.hh"
53 #include "G4Alpha.hh"
54 #include "G4ParticleHPVector.hh"
55 #include "G4NucleiProperties.hh"
57 #include "G4IonTable.hh"
58 #include <set>
59 
61 {
62  theAngular = 0;
63  fCache.Get()->currentMeanEnergy = -2;
64  fCache.Get()->fresh = true;
65  adjustResult = true;
66  if ( getenv( "G4PHP_DO_NOT_ADJUST_FINAL_STATE" ) ) adjustResult = false;
67 
70  theProjectile = projectile;
71 }
72 
73  void G4ParticleHPContAngularPar::Init(std::istream & aDataFile, G4ParticleDefinition* projectile)
74  {
75  adjustResult = true;
76  if ( getenv( "G4PHP_DO_NOT_ADJUST_FINAL_STATE" ) ) adjustResult = false;
77 
78  theProjectile = projectile;
79 
81  /*if( getenv("G4PHPTEST") )*/
82  theEnergy *= eV;
84  for(G4int i=0; i<nEnergies; i++)
85  {
86  G4double sEnergy;
87  aDataFile >> sEnergy;
88  sEnergy*=eV;
89  theAngular[i].SetLabel(sEnergy);
90  theAngular[i].Init(aDataFile, nAngularParameters, 1.);
91  theMinEner = std::min(theMinEner,sEnergy);
92  theMaxEner = std::max(theMaxEner,sEnergy);
93  }
94  }
95 
97  G4ParticleHPContAngularPar::Sample(G4double anEnergy, G4double massCode, G4double /*targetMass*/,
98  G4int angularRep, G4int /*interpolE*/ )
99  {
100  if( getenv("G4PHPTEST") ) G4cout << " G4ParticleHPContAngularPar::Sample " << anEnergy << " " << massCode << " " << angularRep << G4endl; //GDEB
101  if ( fCache.Get() == NULL ) cacheInit();
102  G4ReactionProduct * result = new G4ReactionProduct;
103  G4int Z = static_cast<G4int>(massCode/1000);
104  G4int A = static_cast<G4int>(massCode-1000*Z);
105  if(massCode==0)
106  {
107  result->SetDefinition(G4Gamma::Gamma());
108  }
109  else if(A==0)
110  {
112  if(Z==1) result->SetDefinition(G4Positron::Positron());
113  }
114  else if(A==1)
115  {
117  if(Z==1) result->SetDefinition(G4Proton::Proton());
118  }
119  else if(A==2)
120  {
122  }
123  else if(A==3)
124  {
125  result->SetDefinition(G4Triton::Triton());
126  if(Z==2) result->SetDefinition(G4He3::He3());
127  }
128  else if(A==4)
129  {
130  result->SetDefinition(G4Alpha::Alpha());
131  if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPContAngularPar: Unknown ion case 1");
132  }
133  else
134  {
135  result->SetDefinition(G4IonTable::GetIonTable()->GetIon(Z,A,0));
136  }
137  G4int i(0);
138  G4int it(0);
139  G4double fsEnergy(0);
140  G4double cosTh(0);
141 
142  if( angularRep == 1 )
143  {
144 // 080612 Fix contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #1
145  //if (interpolE == 2)
146 //110609 above was wrong interupition, pointed out by E.Mendoza and D.Cano (CIMAT)
147 //Following are reviesd version written by T.Koi (SLAC)
148  if ( nDiscreteEnergies != 0 )
149  {
150 
151 //1st check remaining_energy
152 // if this is the first set it. (How?)
153  if ( fCache.Get()->fresh == true )
154  {
155  //Discrete Lines, larger energies come first
156  //Continues Emssions, low to high LAST
157  fCache.Get()->remaining_energy = std::max ( theAngular[0].GetLabel() , theAngular[nEnergies-1].GetLabel() );
158  fCache.Get()->fresh = false;
159  }
160 
161  //Cheating for small remaining_energy
162  //TEMPORAL SOLUTION
163  if ( nDiscreteEnergies == nEnergies )
164  {
165  fCache.Get()->remaining_energy = std::max ( fCache.Get()->remaining_energy , theAngular[nDiscreteEnergies-1].GetLabel() ); //Minimum Line
166  }
167  else
168  {
169  //G4double cont_min = theAngular[nDiscreteEnergies].GetLabel();
170  //if ( theAngular[nDiscreteEnergies].GetLabel() == 0.0 ) cont_min = theAngular[nDiscreteEnergies+1].GetLabel();
171  G4double cont_min=0.0;
172  for ( G4int j = nDiscreteEnergies ; j < nEnergies ; j++ )
173  {
174  cont_min = theAngular[j].GetLabel();
175  if ( theAngular[j].GetValue(0) != 0.0 ) break;
176  }
177  fCache.Get()->remaining_energy = std::max ( fCache.Get()->remaining_energy , std::min ( theAngular[nDiscreteEnergies-1].GetLabel() , cont_min ) ); //Minimum Line or grid
178  }
179 //
180  G4double random = G4UniformRand();
181 
182  G4double * running = new G4double[nEnergies+1];
183  running[0] = 0.0;
184 
185  for ( G4int j = 0 ; j < nDiscreteEnergies ; j++ )
186  {
187  G4double delta = 0.0;
188  if ( theAngular[j].GetLabel() <= fCache.Get()->remaining_energy ) delta = theAngular[i].GetValue(0);
189  running[j+1] = running[j] + delta;
190  }
191  G4double tot_prob_DIS = running[ nDiscreteEnergies ];
192 
193  for ( G4int j = nDiscreteEnergies ; j < nEnergies ; j++ )
194  {
195  G4double delta = 0.0;
196  G4double e_low = 0.0;
197  G4double e_high = 0.0;
198  if ( theAngular[j].GetLabel() <= fCache.Get()->remaining_energy ) delta = theAngular[j].GetValue(0);
199 
200  //To calculate Prob. e_low and e_high should be in eV
201  //There are two case
202  //1:theAngular[nDiscreteEnergies].GetLabel() != 0.0
203  // delta should be used between j-1 and j
204  // At j = nDiscreteEnergies (the first) e_low should be set explicitly
205  if ( theAngular[j].GetLabel() != 0 )
206  {
207  if ( j == nDiscreteEnergies ) {
208  e_low = 0.0/eV;
209  } else {
210  e_low = theAngular[j-1].GetLabel()/eV;
211  }
212  e_high = theAngular[j].GetLabel()/eV;
213  }
214  //2:theAngular[nDiscreteEnergies].GetLabel() == 0.0
215  // delta should be used between j and j+1
216  if ( theAngular[j].GetLabel() == 0.0 ) {
217  e_low = theAngular[j].GetLabel()/eV;
218  if ( j != nEnergies-1 ) {
219  e_high = theAngular[j+1].GetLabel()/eV;
220  } else {
221  e_high = theAngular[j].GetLabel()/eV;
222  if ( theAngular[j].GetValue(0) != 0.0 ) {
223  throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPContAngularPar: Unexpected non zero value of theAngular[nEnergies-1].GetValue(0)");
224  }
225  }
226  }
227 
228  running[j+1] = running[j] + ( ( e_high - e_low ) * delta );
229  }
230  G4double tot_prob_CON = running[ nEnergies ] - running[ nDiscreteEnergies ];
231 
232 /*
233  For FPE debugging
234  if (tot_prob_DIS + tot_prob_CON == 0 ) {
235  G4cout << "TKDB tot_prob_DIS + tot_prob_CON " << tot_prob_DIS + tot_prob_CON << G4endl;
236  G4cout << "massCode " << massCode << G4endl;
237  G4cout << "nDiscreteEnergies " << nDiscreteEnergies << " nEnergies " << nEnergies << G4endl;
238  for ( int j = nDiscreteEnergies ; j < nEnergies ; j++ ) {
239  G4cout << j << " " << theAngular[j].GetLabel() << " " << theAngular[j].GetValue(0) << G4endl;
240  }
241  }
242 */
243  // Normalize random
244  random *= (tot_prob_DIS + tot_prob_CON);
245 //2nd Judge Discrete or not This shoudl be relatively close to 1 For safty
246  if ( random <= ( tot_prob_DIS / ( tot_prob_DIS + tot_prob_CON ) ) || nDiscreteEnergies == nEnergies )
247  {
248 // Discrete Emission
249  for ( G4int j = 0 ; j < nDiscreteEnergies ; j++ )
250  {
251  //Here we should use i+1
252  if ( random < running[ j+1 ] )
253  {
254  it = j;
255  break;
256  }
257  }
258  fsEnergy = theAngular[ it ].GetLabel();
259 
260  G4ParticleHPLegendreStore theStore(1);
261  theStore.Init(0,fsEnergy,nAngularParameters);
262  for (G4int j=0;j<nAngularParameters;j++)
263  {
264  theStore.SetCoeff(0,j,theAngular[it].GetValue(j));
265  }
266  // use it to sample.
267  cosTh = theStore.SampleMax(fsEnergy);
268  //Done
269  }
270  else
271  {
272 // Continuous Emission
273  for ( G4int j = nDiscreteEnergies ; j < nEnergies ; j++ )
274  {
275  //Here we should use i
276  if ( random < running[ j ] )
277  {
278  it = j;
279  break;
280  }
281  }
282 
283  G4double x1 = running[it-1];
284  G4double x2 = running[it];
285 
286  G4double y1 = 0.0;
287  if ( it != nDiscreteEnergies )
288  y1 = theAngular[it-1].GetLabel();
289  G4double y2 = theAngular[it].GetLabel();
290 
292  random,x1,x2,y1,y2);
293 
294  G4ParticleHPLegendreStore theStore(2);
295  theStore.Init(0,y1,nAngularParameters);
296  theStore.Init(1,y2,nAngularParameters);
297  theStore.SetManager(theManager);
298  for (G4int j=0;j<nAngularParameters;j++)
299  {
300  G4int itt = it;
301  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
302  if ( it == 0 )
303  {
304  //Safty for unexpected it = 0;
305  //G4cout << "110611 G4ParticleHPContAngularPar::Sample it = 0; invetigation required " << G4endl;
306  itt = it+1;
307  }
308  theStore.SetCoeff(0,j,theAngular[itt-1].GetValue(j));
309  theStore.SetCoeff(1,j,theAngular[itt].GetValue(j));
310  }
311  // use it to sample.
312  cosTh = theStore.SampleMax(fsEnergy);
313 
314  //Done
315  }
316 
317  //TK080711
318  if( adjustResult ) fCache.Get()->remaining_energy -= fsEnergy;
319  //TK080711
320 
321  //080801b
322  delete[] running;
323  //080801b
324  }
325  else
326  {
327  // Only continue, TK will clean up
328 
329  //080714
330  if ( fCache.Get()->fresh == true )
331  {
332  fCache.Get()->remaining_energy = theAngular[ nEnergies-1 ].GetLabel();
333  fCache.Get()->fresh = false;
334  }
335  //080714
336  G4double random = G4UniformRand();
337  G4double * running = new G4double[nEnergies];
338  running[0]=0;
339  G4double weighted = 0;
340  for(i=1; i<nEnergies; i++)
341  {
342 /*
343  if(i!=0)
344  {
345  running[i]=running[i-1];
346  }
347  running[i] += theInt.GetBinIntegral(theManager.GetScheme(i-1),
348  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
349  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
350  weighted += theInt.GetWeightedBinIntegral(theManager.GetScheme(i-1),
351  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
352  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
353 */
354 
355  running[i]=running[i-1];
356  if ( fCache.Get()->remaining_energy >= theAngular[i].GetLabel() )
357  {
358  running[i] += theInt.GetBinIntegral(theManager.GetScheme(i-1),
359  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
360  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
362  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
363  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
364  }
365  }
366  // cash the mean energy in this distribution
367  //080409 TKDB
368  if ( nEnergies == 1 || running[nEnergies-1] == 0 )
369  fCache.Get()->currentMeanEnergy = 0.0;
370  else
371  {
372  fCache.Get()->currentMeanEnergy = weighted/running[nEnergies-1];
373  }
374 
375  //080409 TKDB
376  if ( nEnergies == 1 ) it = 0;
377 
378  //080729
379  if ( running[nEnergies-1] != 0 )
380  {
381  for ( i = 1 ; i < nEnergies ; i++ )
382  {
383  it = i;
384  if ( random < running [ i ] / running [ nEnergies-1 ] ) break;
385  }
386  }
387 
388  //080714
389  if ( running [ nEnergies-1 ] == 0 ) it = 0;
390  //080714
391 
392  if (it<nDiscreteEnergies||it==0)
393  {
394  if(it == 0)
395  {
396  fsEnergy = theAngular[it].GetLabel();
397  G4ParticleHPLegendreStore theStore(1);
398  theStore.Init(0,fsEnergy,nAngularParameters);
399  for(i=0;i<nAngularParameters;i++)
400  {
401  theStore.SetCoeff(0,i,theAngular[it].GetValue(i));
402  }
403  // use it to sample.
404  cosTh = theStore.SampleMax(fsEnergy);
405  }
406  else
407  {
408  G4double e1, e2;
409  e1 = theAngular[it-1].GetLabel();
410  e2 = theAngular[it].GetLabel();
412  random,
413  running[it-1]/running[nEnergies-1],
414  running[it]/running[nEnergies-1],
415  e1, e2);
416  // fill a Legendrestore
417  G4ParticleHPLegendreStore theStore(2);
418  theStore.Init(0,e1,nAngularParameters);
419  theStore.Init(1,e2,nAngularParameters);
420  for(i=0;i<nAngularParameters;i++)
421  {
422  theStore.SetCoeff(0,i,theAngular[it-1].GetValue(i));
423  theStore.SetCoeff(1,i,theAngular[it].GetValue(i));
424  }
425  // use it to sample.
426  theStore.SetManager(theManager);
427  cosTh = theStore.SampleMax(fsEnergy);
428  }
429  }
430  else // continuum contribution
431  {
432  G4double x1 = running[it-1]/running[nEnergies-1];
433  G4double x2 = running[it]/running[nEnergies-1];
434  G4double y1 = theAngular[it-1].GetLabel();
435  G4double y2 = theAngular[it].GetLabel();
437  random,x1,x2,y1,y2);
438  G4ParticleHPLegendreStore theStore(2);
439  theStore.Init(0,y1,nAngularParameters);
440  theStore.Init(1,y2,nAngularParameters);
441  theStore.SetManager(theManager);
442  for(i=0;i<nAngularParameters;i++)
443  {
444  theStore.SetCoeff(0,i,theAngular[it-1].GetValue(i));
445  theStore.SetCoeff(1,i,theAngular[it].GetValue(i));
446  }
447  // use it to sample.
448  cosTh = theStore.SampleMax(fsEnergy);
449  }
450  delete [] running;
451 
452  //080714
453  if( adjustResult ) fCache.Get()->remaining_energy -= fsEnergy;
454  //080714
455  }
456  }
457  else if(angularRep==2)
458  {
459  // first get the energy (already the right for this incoming energy)
460  G4int j;
461  G4double * running = new G4double[nEnergies];
462  running[0]=0;
463  G4double weighted = 0;
464  if( getenv("G4PHPTEST") ) G4cout << " G4ParticleHPContAngularPar::Sample nEnergies " << nEnergies << G4endl;
465  for(j=1; j<nEnergies; j++)
466  {
467  if(j!=0) running[j]=running[j-1];
468  running[j] += theInt.GetBinIntegral(theManager.GetScheme(j-1),
469  theAngular[j-1].GetLabel(), theAngular[j].GetLabel(),
470  theAngular[j-1].GetValue(0), theAngular[j].GetValue(0));
472  theAngular[j-1].GetLabel(), theAngular[j].GetLabel(),
473  theAngular[j-1].GetValue(0), theAngular[j].GetValue(0));
474  if( getenv("G4PHPTEST") ) G4cout << " G4ParticleHPContAngularPar::Sample " << j << " running " << running[j]
475  << " " << theManager.GetScheme(j-1) << " " << theAngular[j-1].GetLabel() << " " << theAngular[j].GetLabel() << " " << theAngular[j-1].GetValue(0) << " " << theAngular[j].GetValue(0) << G4endl; //GDEB
476  }
477  // cash the mean energy in this distribution
478  //080409 TKDB
479  //currentMeanEnergy = weighted/running[nEnergies-1];
480  if ( nEnergies == 1 )
481  fCache.Get()->currentMeanEnergy = 0.0;
482  else
483  fCache.Get()->currentMeanEnergy = weighted/running[nEnergies-1];
484 
485  G4int itt(0);
486  G4double randkal = G4UniformRand();
487  //080409 TKDB
488  //for(i=0; i<nEnergies; i++)
489  for(j=1; j<nEnergies; j++)
490  {
491  itt = j;
492  if(randkal<running[j]/running[nEnergies-1]) break;
493  }
494 
495  // interpolate the secondary energy.
496  G4double x, x1,x2,y1,y2;
497  if(itt==0) itt=1;
498  x = randkal*running[nEnergies-1];
499  x1 = running[itt-1];
500  x2 = running[itt];
501  G4double compoundFraction;
502  // interpolate energy
503  y1 = theAngular[itt-1].GetLabel();
504  y2 = theAngular[itt].GetLabel();
505  fsEnergy = theInt.Interpolate(theManager.GetInverseScheme(itt-1),
506  x, x1,x2,y1,y2);
507  if( getenv("G4PHPTEST") ) G4cout << itt << " G4particleHPContAngularPar fsEnergy " << fsEnergy << " " << theManager.GetInverseScheme(itt-1) << " x " << x << " " << x1 << " " << x2 << " y " << y1 << " " << y2 << G4endl; //GDEB
508  // for theta interpolate the compoundFractions
509  G4double cLow = theAngular[itt-1].GetValue(1);
510  G4double cHigh = theAngular[itt].GetValue(1);
511  compoundFraction = theInt.Interpolate(theManager.GetScheme(itt),
512  fsEnergy, y1, y2, cLow,cHigh);
513  if( getenv("G4PHPTEST") ) G4cout << itt << " G4particleHPContAngularPar compoundFraction " << compoundFraction << " E " << fsEnergy << " " << theManager.GetScheme(itt) << " ener " << fsEnergy << " y " << y1 << " " << y2 << " cLH " << cLow << " " << cHigh << G4endl; //GDEB
514  delete [] running;
515 
516  // get cosTh
517  G4double incidentEnergy = anEnergy;
518  G4double incidentMass = theProjectile->GetPDGMass();
519  G4double productEnergy = fsEnergy;
520  G4double productMass = result->GetMass();
521  G4int targetZ = G4int(theTargetCode/1000);
522  G4int targetA = G4int(theTargetCode-1000*targetZ);
523  // To correspond to natural composition (-nat-) data files.
524  if ( targetA == 0 )
525  targetA = G4int ( theTarget->GetMass()/amu_c2 + 0.5 );
526  G4double targetMass = theTarget->GetMass();
527  G4int residualA = targetA+1-A;
528  G4int residualZ = targetZ-Z;
529  G4double residualMass = residualZ*G4Proton::Proton()->GetPDGMass();
530  residualMass +=(residualA-residualZ)*theProjectile->GetPDGMass();
531  residualMass -= G4NucleiProperties::GetBindingEnergy( residualA , residualZ );
532  G4ParticleHPKallbachMannSyst theKallbach(compoundFraction,
533  incidentEnergy, incidentMass,
534  productEnergy, productMass,
535  residualMass, residualA, residualZ,
536  targetMass, targetA, targetZ);
537  cosTh = theKallbach.Sample(anEnergy);
538  if( getenv("G4PHPTEST") ) G4cout << " G4ParticleHPKallbachMannSyst::Sample resulttest " << cosTh << G4endl; //GDEB
539  }
540  else if(angularRep>10&&angularRep<16)
541  {
542  G4double random = G4UniformRand();
543  G4double * running = new G4double[nEnergies];
544  running[0]=0;
545  G4double weighted = 0;
546  for(i=1; i<nEnergies; i++)
547  {
548  if(i!=0) running[i]=running[i-1];
549  running[i] += theInt.GetBinIntegral(theManager.GetScheme(i-1),
550  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
551  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
553  theAngular[i-1].GetLabel(), theAngular[i].GetLabel(),
554  theAngular[i-1].GetValue(0), theAngular[i].GetValue(0));
555  }
556  // cash the mean energy in this distribution
557  //currentMeanEnergy = weighted/running[nEnergies-1];
558  if ( nEnergies == 1 )
559  fCache.Get()->currentMeanEnergy = 0.0;
560  else
561  fCache.Get()->currentMeanEnergy = weighted/running[nEnergies-1];
562 
563  //080409 TKDB
564  if ( nEnergies == 1 ) it = 0;
565  //for(i=0; i<nEnergies; i++)
566  for(i=1; i<nEnergies; i++)
567  {
568  it = i;
569  if(random<running[i]/running[nEnergies-1]) break;
570  }
571  if(it<nDiscreteEnergies||it==0)
572  {
573  if(it==0)
574  {
575  fsEnergy = theAngular[0].GetLabel();
576  G4ParticleHPVector theStore;
577  G4int aCounter = 0;
578  for(G4int j=1; j<nAngularParameters; j+=2)
579  {
580  theStore.SetX(aCounter, theAngular[0].GetValue(j));
581  theStore.SetY(aCounter, theAngular[0].GetValue(j+1));
582  aCounter++;
583  }
585  aMan.Init(angularRep-10, nAngularParameters-1);
586  theStore.SetInterpolationManager(aMan);
587  cosTh = theStore.Sample();
588  }
589  else
590  {
591  fsEnergy = theAngular[it].GetLabel();
592  G4ParticleHPVector theStore;
594  aMan.Init(angularRep-10, nAngularParameters-1);
595  theStore.SetInterpolationManager(aMan); // Store interpolates f(costh)
597  G4int aCounter = 0;
598  for(G4int j=1; j<nAngularParameters; j+=2)
599  {
600  theStore.SetX(aCounter, theAngular[it].GetValue(j));
601  theStore.SetY(aCounter, theInt.Interpolate(currentScheme,
602  random,
603  running[it-1]/running[nEnergies-1],
604  running[it]/running[nEnergies-1],
605  theAngular[it-1].GetValue(j+1),
606  theAngular[it].GetValue(j+1)));
607  aCounter++;
608  }
609  cosTh = theStore.Sample();
610  }
611  }
612  else
613  {
614  G4double x1 = running[it-1]/running[nEnergies-1];
615  G4double x2 = running[it]/running[nEnergies-1];
616  G4double y1 = theAngular[it-1].GetLabel();
617  G4double y2 = theAngular[it].GetLabel();
619  random,x1,x2,y1,y2);
620  G4ParticleHPVector theBuff1;
621  G4ParticleHPVector theBuff2;
623  aMan.Init(angularRep-10, nAngularParameters-1);
624 // theBuff1.SetInterpolationManager(aMan); // Store interpolates f(costh)
625 // theBuff2.SetInterpolationManager(aMan); // Store interpolates f(costh)
626 // Bug Report #1366 from L. Russell
627  //for(i=0; i<nAngularParameters; i++) // i=1 ist wichtig!
628  //{
629  // theBuff1.SetX(i, theAngular[it-1].GetValue(i));
630  // theBuff1.SetY(i, theAngular[it-1].GetValue(i+1));
631  // theBuff2.SetX(i, theAngular[it].GetValue(i));
632  // theBuff2.SetY(i, theAngular[it].GetValue(i+1));
633  // i++;
634  //}
635  {
636  G4int j;
637  for(i=0,j=1; i<nAngularParameters; i++,j+=2)
638  {
639  theBuff1.SetX(i, theAngular[it-1].GetValue(j));
640  theBuff1.SetY(i, theAngular[it-1].GetValue(j+1));
641  theBuff2.SetX(i, theAngular[it].GetValue(j));
642  theBuff2.SetY(i, theAngular[it].GetValue(j+1));
643  }
644  }
645  G4ParticleHPVector theStore;
646  theStore.SetInterpolationManager(aMan); // Store interpolates f(costh)
647  x1 = y1;
648  x2 = y2;
649  G4double x, y;
650  //for(i=0;i<theBuff1.GetVectorLength(); i++);
651  for(i=0;i<theBuff1.GetVectorLength(); i++)
652  {
653  x = theBuff1.GetX(i); // costh binning identical
654  y1 = theBuff1.GetY(i);
655  y2 = theBuff2.GetY(i);
657  fsEnergy, theAngular[it-1].GetLabel(),
658  theAngular[it].GetLabel(), y1, y2);
659  theStore.SetX(i, x);
660  theStore.SetY(i, y);
661  }
662  cosTh = theStore.Sample();
663  }
664  delete [] running;
665  }
666  else
667  {
668  throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPContAngularPar::Sample: Unknown angular representation");
669  }
670  result->SetKineticEnergy(fsEnergy);
671  G4double phi = twopi*G4UniformRand();
672  G4double theta = std::acos(cosTh);
673  G4double sinth = std::sin(theta);
674  G4double mtot = result->GetTotalMomentum();
675  G4ThreeVector tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) );
676  result->SetMomentum(tempVector);
677 // return the result.
678  return result;
679  }
680 
681 
682 #define MERGE_NEW
683 
685 {
686 
687  //----- Discrete energies: store own energies in a map for faster searching
688  G4int ie;
689  for(ie=0; ie<nDiscreteEnergies; ie++) {
691  }
692  if( !angParPrev ) return;
693 
694  //----- Discrete energies: use energies that appear in one or another
695  for(ie=0; ie<nDiscreteEnergies; ie++) {
696  theDiscreteEnergies.insert(theAngular[ie].GetLabel());
697  }
698  G4int nDiscreteEnergiesPrev = angParPrev->GetNDiscreteEnergies();
699  for(ie=0; ie<nDiscreteEnergiesPrev; ie++) {
700  theDiscreteEnergies.insert(angParPrev->theAngular[ie].GetLabel());
701  }
702 
703  //--- Get the values for which interpolation will be done : all energies of this and previous ContAngularPar
704  for(ie=nDiscreteEnergies; ie<nEnergies; ie++) {
705  G4double ener = theAngular[ie].GetLabel();
706  G4double enerT = (ener-theMinEner)/(theMaxEner-theMinEner);
707  theEnergiesTransformed.insert(enerT);
708  //- if( getenv("G4PHPTEST2") ) G4cout <<this << " G4ParticleHPContAngularPar::PrepareTableInterpolation theEnergiesTransformed1 " << enerT << G4endl; //GDEB
709  }
710  G4int nEnergiesPrev = angParPrev->GetNEnergies();
711  G4double minEnerPrev = angParPrev->GetMinEner();
712  G4double maxEnerPrev = angParPrev->GetMaxEner();
713  for(ie=nDiscreteEnergiesPrev; ie<nEnergiesPrev; ie++) {
714  G4double ener = angParPrev->theAngular[ie].GetLabel();
715  G4double enerT = (ener-minEnerPrev)/(maxEnerPrev-minEnerPrev);
716  theEnergiesTransformed.insert(enerT);
717  //- if( getenv("G4PHPTEST2") ) G4cout << this << " G4ParticleHPContAngularPar::PrepareTableInterpolation theEnergiesTransformed2 " << enerT << G4endl; //GDEB
718  }
719  // add the maximum energy
720  theEnergiesTransformed.insert(1.);
721 
722 }
723 
725  G4ParticleHPContAngularPar & angpar1,
726  G4ParticleHPContAngularPar & angpar2)
727 {
728  G4int ie,ie1,ie2, ie1Prev, ie2Prev;
730  theManager = angpar1.theManager;
731  theEnergy = anEnergy;
732 
734  std::set<G4double>::const_iterator itede;
735  std::map<G4double,G4int> discEnerOwn1 = angpar1.GetDiscreteEnergiesOwn();
736  std::map<G4double,G4int> discEnerOwn2 = angpar2.GetDiscreteEnergiesOwn();
737  std::map<G4double,G4int>::const_iterator itedeo;
738  ie = 0;
739  for( itede = theDiscreteEnergies.begin(); itede != theDiscreteEnergies.end(); itede++, ie++ ) {
740  G4double discEner = *itede;
741  itedeo = discEnerOwn1.find(discEner);
742  if( itedeo == discEnerOwn1.end() ) {
743  ie1 = 0;
744  } else {
745  ie1 = -1;
746  }
747  itedeo = discEnerOwn2.find(discEner);
748  if( itedeo == discEnerOwn2.end() ) {
749  ie2 = 0;
750  } else {
751  ie2 = -1;
752  }
753 
754  theAngular[ie].SetLabel(discEner);
755  G4double val1, val2;
756  for(G4int ip=0; ip<nAngularParameters; ip++) {
757  if( ie1 != -1 ) {
758  val1 = angpar1.theAngular[ie1].GetValue(ip);
759  } else {
760  val1 = 0.;
761  }
762  if( ie2 != -1 ) {
763  val2 = angpar2.theAngular[ie2].GetValue(ip);
764  } else {
765  val2 = 0.;
766  }
767 
768  G4double value = theInt.Interpolate(aScheme, anEnergy,
769  angpar1.theEnergy, angpar2.theEnergy,
770  val1,
771  val2);
772  if( getenv("G4PHPTEST2") ) G4cout << ie << " " << ip << " G4ParticleHPContAngularPar::Merge DiscreteEnergies val1 " << val1 << " val2 " << val2 << " value " << value << G4endl; //GDEB
773 
774  theAngular[ie].SetValue(ip, value);
775  }
776  }
777 
778  if(theAngular != 0) delete [] theAngular;
781 
782  //---- Get minimum and maximum energy interpolating
783  theMinEner = angpar1.GetMinEner() + (theEnergy-angpar1.GetEnergy()) * (angpar2.GetMinEner()-angpar1.GetMinEner())/(angpar2.GetEnergy()-angpar1.GetEnergy());
784  theMaxEner = angpar1.GetMaxEner() + (theEnergy-angpar1.GetEnergy()) * (angpar2.GetMaxEner()-angpar1.GetMaxEner())/(angpar2.GetEnergy()-angpar1.GetEnergy());
785 
786  if( getenv("G4PHPTEST2") ) G4cout << " G4ParticleHPContAngularPar::Merge E " << anEnergy << " minmax " << theMinEner << " " << theMaxEner << G4endl; //GDEB
787 
788  //--- Loop to energies of new set
789  std::set<G4double> energiesTransformed = angpar2.GetEnergiesTransformed();
790  std::set<G4double>::const_iterator iteet = energiesTransformed.begin();
791  G4int nEnergies1 = angpar1.GetNEnergies();
792  G4int nDiscreteEnergies1 = angpar1.GetNDiscreteEnergies();
793  G4double minEner1 = angpar1.GetMinEner();
794  G4double maxEner1 = angpar1.GetMaxEner();
795  G4int nEnergies2 = angpar2.GetNEnergies();
796  G4int nDiscreteEnergies2 = angpar2.GetNDiscreteEnergies();
797  G4double minEner2 = angpar2.GetMinEner();
798  G4double maxEner2 = angpar2.GetMaxEner();
799  for(ie=nDiscreteEnergies; ie<nEnergies; ie++,iteet++) {
800  G4double eT = (*iteet);
801 
802  //--- Use eT1 = eT: Get energy and parameters of angpar1 for this eT
803  G4double e1 = (maxEner1-minEner1) * eT + minEner1;
804  //----- Get parameter value corresponding to this e1
805  for(ie1=nDiscreteEnergies1; ie1<nEnergies1; ie1++) {
806  if( (angpar1.theAngular[ie1].GetLabel() - e1) > 1.E-10*e1 ) break;
807  }
808  ie1Prev = ie1 - 1;
809  if( ie1 == 0 ) ie1Prev++;
810  if( ie1 == nEnergies1 ) {
811  ie1--;
812  ie1Prev = ie1;
813  }
814  //--- Use eT2 = eT: Get energy and parameters of angpar2 for this eT
815  G4double e2 = (maxEner2-minEner2) * eT + minEner2;
816  //----- Get parameter value corresponding to this e2
817  for(ie2=nDiscreteEnergies2; ie2<nEnergies2; ie2++) {
818  // G4cout << " GET IE2 " << ie2 << " - " << angpar2.theAngular[ie2].GetLabel() - e2 << " " << angpar2.theAngular[ie2].GetLabel() << " " << e2 <<G4endl;
819  if( (angpar2.theAngular[ie2].GetLabel() - e2) > 1.E-10*e2 ) break;
820  }
821  ie2Prev = ie2 - 1;
822  if( ie2 == 0 ) ie2Prev++;
823  if( ie2 == nEnergies2 ) {
824  ie2--;
825  ie2Prev = ie2;
826  }
827 
828  //---- Energy corresponding to energy transformed
830  if( getenv("G4PHPTEST2") ) G4cout << ie << " " << ie1 << " " << ie2 << " G4ParticleHPContAngularPar::loop eT " << eT << " -> eN " << eN << " e1 " << e1 << " e2 " << e2 << G4endl; //GDEB
831 
832  theAngular[ie].SetLabel(eN);
833 
834  for(G4int ip=0; ip<nAngularParameters; ip++) {
836  e1,
837  angpar1.theAngular[ie1Prev].GetLabel(),
838  angpar1.theAngular[ie1].GetLabel(),
839  angpar1.theAngular[ie1Prev].GetValue(ip),
840  angpar1.theAngular[ie1].GetValue(ip)) * (maxEner1-minEner1);
842  e2,
843  angpar2.theAngular[ie2Prev].GetLabel(),
844  angpar2.theAngular[ie2].GetLabel(),
845  angpar2.theAngular[ie2Prev].GetValue(ip),
846  angpar2.theAngular[ie2].GetValue(ip)) * (maxEner2-minEner2);
847 
848  G4double value = theInt.Interpolate(aScheme, anEnergy,
849  angpar1.theEnergy, angpar2.theEnergy,
850  val1,
851  val2);
852  value /= (theMaxEner-theMinEner);
853  if( getenv("G4PHPTEST2") ) G4cout << ie << " " << ip << " G4ParticleHPContAngularPar::Merge val1 " << val1 << " val2 " << val2 << " value " << value << G4endl; //GDEB
854  //- val1 = angpar1.theAngular[ie1-1].GetValue(ip) * (maxEner1-minEner1);
855  //- val2 = angpar2.theAngular[ie2-1].GetValue(ip) * (maxEner2-minEner2);
856  //- if( getenv("G4PHPTEST2") ) G4cout << ie << " " << ip << " G4ParticleHPContAngularPar::MergeOLD val1 " << val1 << " val2 " << val2 << " value " << value << G4endl; //GDEB
857 
858  theAngular[ie].SetValue(ip, value);
859  }
860  }
861 
862  if( getenv("G4PHPTEST2") ) {
863  G4cout << " G4ParticleHPContAngularPar::Merge ANGPAR1 " << G4endl; //GDEB
864  angpar1.Dump();
865  G4cout << " G4ParticleHPContAngularPar::Merge ANGPAR2 " << G4endl;
866  angpar2.Dump();
867  G4cout << " G4ParticleHPContAngularPar::Merge ANGPARNEW " << G4endl;
868  Dump();
869  }
870 }
871 
873 {
874  G4cout << theEnergy << " " << nEnergies << " " << nDiscreteEnergies << " " << nAngularParameters << G4endl;
875 
876  for(G4int ii=0; ii<nEnergies; ii++) {
877  theAngular[ii].Dump();
878  }
879 
880 }
std::map< G4double, G4int > GetDiscreteEnergiesOwn() const
void Init(G4int i, G4double e, G4int n)
G4int GetVectorLength() const
G4double GetTotalMomentum() const
CLHEP::Hep3Vector G4ThreeVector
void Init(std::istream &aDataFile, G4int nPar, G4double unit=1.)
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass, G4int angularRep, G4int interpol)
void BuildByInterpolation(G4double anEnergy, G4InterpolationScheme aScheme, G4ParticleHPContAngularPar &store1, G4ParticleHPContAngularPar &store2)
G4double Interpolate2(G4InterpolationScheme aScheme, G4double x, G4double x1, G4double x2, G4double y1, G4double y2) const
void Init(G4int aScheme, G4int aRange)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double SampleMax(G4double energy)
static const G4double e2
void SetInterpolationManager(const G4InterpolationManager &aManager)
int G4int
Definition: G4Types.hh:78
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
std::map< G4double, G4int > theDiscreteEnergiesOwn
std::set< G4double > GetEnergiesTransformed() const
void SetY(G4int i, G4double x)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
void Init(std::istream &aDataFile, G4ParticleDefinition *projectile)
void SetValue(G4int i, G4double y)
G4InterpolationScheme GetScheme(G4int index) const
static const double twopi
Definition: G4SIunits.hh:75
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetX(G4int i) const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetY(G4double x)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetWeightedBinIntegral(const G4InterpolationScheme &aScheme, const G4double x1, const G4double x2, const G4double y1, const G4double y2)
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double Interpolate(G4InterpolationScheme aScheme, G4double x, G4double x1, G4double x2, G4double y1, G4double y2) const
void SetLabel(G4double aLabel)
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static const G4double e1
void SetManager(G4InterpolationManager &aManager)
static const double eV
Definition: G4SIunits.hh:212
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4InterpolationScheme
G4double GetPDGMass() const
static G4double GetBindingEnergy(const G4int A, const G4int Z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetCoeff(G4int i, G4int l, G4double coeff)
const G4double x[NPOINTSGL]
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4InterpolationScheme GetInverseScheme(G4int index)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetValue(G4int i)
#define G4endl
Definition: G4ios.hh:61
G4double GetBinIntegral(const G4InterpolationScheme &aScheme, const G4double x1, const G4double x2, const G4double y1, const G4double y2)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
double G4double
Definition: G4Types.hh:76
void PrepareTableInterpolation(const G4ParticleHPContAngularPar *angularPrev)
static G4He3 * He3()
Definition: G4He3.cc:94
void SetX(G4int i, G4double e)
#define DBL_MAX
Definition: templates.hh:83
G4double GetMass() const