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