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