Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPThermalScattering.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// Thermal Neutron Scattering
27// Koi, Tatsumi (SLAC/SCCS)
28//
29// Class Description:
30//
31// Final State Generators for a high precision (based on evaluated data
32// libraries) description of themal neutron scattering below 4 eV;
33// Based on Thermal neutron scattering files
34// from the evaluated nuclear data files ENDF/B-VI, Release2
35// To be used in your physics list in case you need this physics.
36// In this case you want to register an object of this class with
37// the corresponding process.
38
39
40// 070625 Fix memory leaking at destructor by T. Koi
41// 081201 Fix memory leaking at destructor by T. Koi
42// 100729 Add model name in constructor Problem #1116
43// P. Arce, June-2014 Conversion neutron_hp to particle_hp
44//
50#include "G4SystemOfUnits.hh"
51#include "G4Neutron.hh"
52#include "G4ElementTable.hh"
53#include "G4MaterialTable.hh"
54#include "G4Threading.hh"
55
57 :G4HadronicInteraction("NeutronHPThermalScattering")
58,coherentFSs(nullptr)
59,incoherentFSs(nullptr)
60,inelasticFSs(nullptr)
61{
62 theHPElastic = new G4ParticleHPElastic();
63
64 SetMinEnergy( 0.*eV );
65 SetMaxEnergy( 4*eV );
66 theXSection = new G4ParticleHPThermalScatteringData();
67
68 nMaterial = 0;
69 nElement = 0;
70}
71
72
74{
75 delete theHPElastic;
76}
77
78
79void G4ParticleHPThermalScattering::clearCurrentFSData()
80{
81 if ( incoherentFSs != nullptr )
82 {
83 for (auto it=incoherentFSs->cbegin(); it!=incoherentFSs->cend(); ++it)
84 {
85 for (auto itt=it->second->cbegin(); itt!= it->second->cend(); ++itt)
86 {
87 for (auto ittt=itt->second->cbegin(); ittt!=itt->second->cend(); ++ittt)
88 {
89 delete *ittt;
90 }
91 delete itt->second;
92 }
93 delete it->second;
94 }
95 }
96
97 if ( coherentFSs != nullptr )
98 {
99 for (auto it=coherentFSs->cbegin(); it!=coherentFSs->cend(); ++it)
100 {
101 for (auto itt=it->second->cbegin(); itt!=it->second->cend(); ++itt)
102 {
103 for (auto ittt=itt->second->cbegin(); ittt!=itt->second->cend(); ++ittt)
104 {
105 delete *ittt;
106 }
107 delete itt->second;
108 }
109 delete it->second;
110 }
111 }
112
113 if ( inelasticFSs != nullptr )
114 {
115 for (auto it=inelasticFSs->cbegin(); it!=inelasticFSs->cend(); ++it)
116 {
117 for (auto itt=it->second->cbegin(); itt!=it->second->cend(); ++itt)
118 {
119 for (auto ittt=itt->second->cbegin(); ittt!=itt->second->cend(); ++ittt)
120 {
121 for (auto it4=(*ittt)->vE_isoAngle.cbegin(); it4!=(*ittt)->vE_isoAngle.cend(); ++it4)
122 {
123 delete *it4;
124 }
125 delete *ittt;
126 }
127 delete itt->second;
128 }
129 delete it->second;
130 }
131 }
132
133 incoherentFSs = nullptr;
134 coherentFSs = nullptr;
135 inelasticFSs = nullptr;
136}
137
138
140{
141 buildPhysicsTable();
142 theHPElastic->BuildPhysicsTable( particle );
143}
144
145
146std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >*
147G4ParticleHPThermalScattering::readACoherentFSDATA( G4String name )
148{
149 auto aCoherentFSDATA = new std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >;
150
151 std::istringstream theChannel(std::ios::in);
153
154 std::vector< G4double > vBraggE;
155
156 G4int dummy;
157 while ( theChannel >> dummy ) // MF // Loop checking, 11.05.2015, T. Koi
158 {
159 theChannel >> dummy; // MT
160 G4double temp;
161 theChannel >> temp;
162 std::vector < std::pair< G4double , G4double >* >*
163 anBragE_P = new std::vector < std::pair< G4double , G4double >* >;
164
165 G4int n;
166 theChannel >> n;
167 for ( G4int i = 0 ; i < n ; ++i )
168 {
169 G4double Ei;
170 G4double Pi;
171 if ( aCoherentFSDATA->size() == 0 )
172 {
173 theChannel >> Ei;
174 vBraggE.push_back( Ei );
175 }
176 else
177 {
178 Ei = vBraggE[ i ];
179 }
180 theChannel >> Pi;
181 anBragE_P->push_back ( new std::pair < G4double , G4double > ( Ei , Pi ) );
182 }
183 aCoherentFSDATA->insert ( std::pair < G4double , std::vector < std::pair< G4double , G4double >* >* > ( temp , anBragE_P ) );
184 }
185 return aCoherentFSDATA;
186}
187
188
189std::map < G4double , std::vector < E_P_E_isoAng* >* >*
190G4ParticleHPThermalScattering::readAnInelasticFSDATA ( G4String name )
191{
192 auto anT_E_P_E_isoAng = new std::map < G4double , std::vector < E_P_E_isoAng* >* >;
193
194 std::istringstream theChannel(std::ios::in);
196
197 G4int dummy;
198 while ( theChannel >> dummy ) // MF // Loop checking, 11.05.2015, T. Koi
199 {
200 theChannel >> dummy; // MT
201 G4double temp;
202 theChannel >> temp;
203 std::vector < E_P_E_isoAng* >* vE_P_E_isoAng = new std::vector < E_P_E_isoAng* >;
204 G4int n;
205 theChannel >> n;
206 for ( G4int i = 0 ; i < n ; ++i )
207 {
208 vE_P_E_isoAng->push_back ( readAnE_P_E_isoAng ( &theChannel ) );
209 }
210 anT_E_P_E_isoAng->insert ( std::pair < G4double , std::vector < E_P_E_isoAng* >* > ( temp , vE_P_E_isoAng ) );
211 }
212
213 return anT_E_P_E_isoAng;
214}
215
216
218G4ParticleHPThermalScattering::readAnE_P_E_isoAng( std::istream* file ) // for inelastic
219{
220 E_P_E_isoAng* aData = new E_P_E_isoAng;
221
222 G4double dummy;
224 G4int nep , nl;
225 *file >> dummy;
226 *file >> energy;
227 aData->energy = energy*eV;
228 *file >> dummy;
229 *file >> dummy;
230 *file >> nep;
231 *file >> nl;
232 aData->n = nep/nl;
233 for ( G4int i = 0 ; i < aData->n ; ++i )
234 {
235 G4double prob;
236 E_isoAng* anE_isoAng = new E_isoAng;
237 aData->vE_isoAngle.push_back( anE_isoAng );
238 *file >> energy;
239 anE_isoAng->energy = energy*eV;
240 anE_isoAng->n = nl - 2;
241 anE_isoAng->isoAngle.resize( anE_isoAng->n );
242 *file >> prob;
243 aData->prob.push_back( prob );
244 //G4cout << "G4ParticleHPThermalScattering inelastic " << energy/eV << " " << i << " " << prob << " " << aData->prob[ i ] << G4endl;
245 for ( G4int j = 0 ; j < anE_isoAng->n ; ++j )
246 {
247 G4double x;
248 *file >> x;
249 anE_isoAng->isoAngle[j] = x ;
250 }
251 }
252
253 // Calcuate sum_of_provXdEs
254 G4double total = 0;
255 aData->secondary_energy_cdf.push_back(0.);
256 for ( G4int i = 0 ; i < aData->n - 1 ; ++i )
257 {
258 G4double E_L = aData->vE_isoAngle[i]->energy/eV;
259 G4double E_H = aData->vE_isoAngle[i+1]->energy/eV;
260 G4double dE = E_H - E_L;
261 G4double pdf = (aData->prob[i] + aData->prob[i+1] )/2. * dE;
262 total += ( pdf );
263 aData->secondary_energy_cdf.push_back( total );
264 aData->secondary_energy_pdf.push_back( pdf );
265 aData->secondary_energy_value.push_back( E_L );
266 }
267
268 aData->sum_of_probXdEs = total;
269
270 // Normalize CDF
272 for ( G4int i = 0; i < aData->secondary_energy_cdf_size; ++i )
273 {
274 aData->secondary_energy_cdf[i] /= total;
275 }
276
277 return aData;
278}
279
280
281std::map < G4double , std::vector < E_isoAng* >* >*
282G4ParticleHPThermalScattering::readAnIncoherentFSDATA ( G4String name )
283{
284 auto T_E = new std::map < G4double , std::vector < E_isoAng* >* >;
285
286 //std::ifstream theChannel( name.c_str() );
287 std::istringstream theChannel(std::ios::in);
289
290 G4int dummy;
291 while ( theChannel >> dummy ) // MF // Loop checking, 11.05.2015, T. Koi
292 {
293 theChannel >> dummy; // MT
294 G4double temp;
295 theChannel >> temp;
296 std::vector < E_isoAng* >* vE_isoAng = new std::vector < E_isoAng* >;
297 G4int n;
298 theChannel >> n;
299 for ( G4int i = 0 ; i < n ; i++ )
300 vE_isoAng->push_back ( readAnE_isoAng( &theChannel ) );
301 T_E->insert ( std::pair < G4double , std::vector < E_isoAng* >* > ( temp , vE_isoAng ) );
302 }
303 //theChannel.close();
304
305 return T_E;
306}
307
308
309E_isoAng* G4ParticleHPThermalScattering::readAnE_isoAng( std::istream* file )
310{
311 E_isoAng* aData = new E_isoAng;
312
313 G4double dummy;
315 G4int n;
316 *file >> dummy;
317 *file >> energy;
318 *file >> dummy;
319 *file >> dummy;
320 *file >> n;
321 *file >> dummy;
322 aData->energy = energy*eV;
323 aData->n = n-2;
324 aData->isoAngle.resize( n );
325
326 *file >> dummy;
327 *file >> dummy;
328 for ( G4int i = 0 ; i < aData->n ; i++ )
329 *file >> aData->isoAngle[i];
330
331 return aData;
332}
333
334
336{
337
338// Select Element > Reaction >
339
340 const G4Material * theMaterial = aTrack.GetMaterial();
341 G4double aTemp = theMaterial->GetTemperature();
342 G4int n = (G4int)theMaterial->GetNumberOfElements();
343
344 G4bool findThermalElement = false;
345 G4int ielement;
346 const G4Element* theElement = nullptr;
347 for ( G4int i = 0; i < n ; ++i )
348 {
349 theElement = theMaterial->GetElement(i);
350 // Select target element
351 if ( aNucleus.GetZ_asInt() == (G4int)(theElement->GetZ() + 0.5 ) )
352 {
353 //Check Applicability of Thermal Scattering
354 if ( getTS_ID( nullptr , theElement ) != -1 )
355 {
356 ielement = getTS_ID( nullptr , theElement );
357 findThermalElement = true;
358 break;
359 }
360 else if ( getTS_ID( theMaterial , theElement ) != -1 )
361 {
362 ielement = getTS_ID( theMaterial , theElement );
363 findThermalElement = true;
364 break;
365 }
366 }
367 }
368
369 if ( findThermalElement == true )
370 {
371
372 // Select Reaction (Inelastic, coherent, incoherent)
373 const G4ParticleDefinition* pd = aTrack.GetDefinition();
374 G4DynamicParticle* dp = new G4DynamicParticle ( pd , aTrack.Get4Momentum() );
375 G4double total = theXSection->GetCrossSection( dp , theElement , theMaterial );
376 G4double inelastic = theXSection->GetInelasticCrossSection( dp , theElement , theMaterial );
377
378 G4double random = G4UniformRand();
379 if ( random <= inelastic/total )
380 {
381 // Inelastic
382
383 std::vector<G4double> v_temp;
384 v_temp.clear();
385 for (auto it = inelasticFSs->find( ielement )->second->cbegin();
386 it != inelasticFSs->find( ielement )->second->cend() ; ++it )
387 {
388 v_temp.push_back( it->first );
389 }
390
391 std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
392 //
393 // For T_L aNEP_EPM_TL and T_H aNEP_EPM_TH
394 //
395 std::vector< E_P_E_isoAng* >* vNEP_EPM_TL = nullptr;
396 std::vector< E_P_E_isoAng* >* vNEP_EPM_TH = nullptr;
397
398 if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
399 {
400 vNEP_EPM_TL = inelasticFSs->find( ielement )->second->find ( tempLH.first/kelvin )->second;
401 vNEP_EPM_TH = inelasticFSs->find( ielement )->second->find ( tempLH.second/kelvin )->second;
402 }
403 else if ( tempLH.first == 0.0 )
404 {
405 auto itm = inelasticFSs->find( ielement )->second->cbegin();
406 vNEP_EPM_TL = itm->second;
407 ++itm;
408 vNEP_EPM_TH = itm->second;
409 tempLH.first = tempLH.second;
410 tempLH.second = itm->first;
411 }
412 else if ( tempLH.second == 0.0 )
413 {
414 auto itm = inelasticFSs->find( ielement )->second->cend();
415 --itm;
416 vNEP_EPM_TH = itm->second;
417 --itm;
418 vNEP_EPM_TL = itm->second;
419 tempLH.second = tempLH.first;
420 tempLH.first = itm->first;
421 }
422
423 G4double sE=0., mu=1.0;
424
425 // New Geant4 method - Stochastic temperature interpolation of the final state
426 // (continuous temperature interpolation was used previously)
427 std::pair< G4double , G4double > secondaryParam;
428 G4double rand_temp = G4UniformRand();
429 if ( rand_temp < (aTemp-tempLH.first)/(tempLH.second - tempLH.first) )
430 secondaryParam = sample_inelastic_E_mu( aTrack.GetKineticEnergy() , vNEP_EPM_TH );
431 else
432 secondaryParam = sample_inelastic_E_mu( aTrack.GetKineticEnergy() , vNEP_EPM_TL );
433
434 sE = secondaryParam.first;
435 mu = secondaryParam.second;
436
437 //set
439 G4double phi = CLHEP::twopi*G4UniformRand();
440 G4double sint= std::sqrt ( 1 - mu*mu );
441 theParticleChange.SetMomentumChange( sint*std::cos(phi), sint*std::sin(phi), mu );
442 }
443 else if ( random <= ( inelastic + theXSection->GetCoherentCrossSection( dp , theElement , theMaterial ) ) / total )
444 {
445 // Coherent Elastic
446
447 G4double E = aTrack.GetKineticEnergy();
448
449 // T_L and T_H
450 std::vector<G4double> v_temp;
451 v_temp.clear();
452 for (auto it = coherentFSs->find(ielement)->second->cbegin();
453 it != coherentFSs->find(ielement)->second->cend(); ++it)
454 {
455 v_temp.push_back( it->first );
456 }
457
458 // T_L T_H
459 std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
460 //
461 //
462 // For T_L anEPM_TL and T_H anEPM_TH
463 //
464 std::vector< std::pair< G4double , G4double >* >* pvE_p_TL = nullptr;
465 std::vector< std::pair< G4double , G4double >* >* pvE_p_TH = nullptr;
466
467 if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
468 {
469 pvE_p_TL = coherentFSs->find( ielement )->second->find ( tempLH.first/kelvin )->second;
470 pvE_p_TH = coherentFSs->find( ielement )->second->find ( tempLH.first/kelvin )->second;
471 }
472 else if ( tempLH.first == 0.0 )
473 {
474 pvE_p_TL = coherentFSs->find( ielement )->second->find ( v_temp[ 0 ] )->second;
475 pvE_p_TH = coherentFSs->find( ielement )->second->find ( v_temp[ 1 ] )->second;
476 tempLH.first = tempLH.second;
477 tempLH.second = v_temp[ 1 ];
478 }
479 else if ( tempLH.second == 0.0 )
480 {
481 pvE_p_TH = coherentFSs->find( ielement )->second->find ( v_temp.back() )->second;
482 auto itv = v_temp.cend();
483 --itv;
484 --itv;
485 pvE_p_TL = coherentFSs->find( ielement )->second->find ( *itv )->second;
486 tempLH.second = tempLH.first;
487 tempLH.first = *itv;
488 }
489 else
490 {
491 // tempLH.first == 0.0 && tempLH.second
492 throw G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Unexpected temperature values in data");
493 }
494
495 std::vector< G4double > vE_T;
496 std::vector< G4double > vp_T;
497
498 G4int n1 = (G4int)pvE_p_TL->size();
499
500 // New Geant4 method - Stochastic interpolation of the final state
501 std::vector< std::pair< G4double , G4double >* >* pvE_p_T_sampled;
502 G4double rand_temp = G4UniformRand();
503 if ( rand_temp < (aTemp-tempLH.first)/(tempLH.second - tempLH.first) )
504 pvE_p_T_sampled = pvE_p_TH;
505 else
506 pvE_p_T_sampled = pvE_p_TL;
507
508 //171005 fix bug, contribution from H.N. TRAN@CEA
509 for ( G4int i=0 ; i < n1 ; ++i )
510 {
511 vE_T.push_back ( (*pvE_p_T_sampled)[i]->first );
512 vp_T.push_back ( (*pvE_p_T_sampled)[i]->second );
513 }
514
515 G4int j = 0;
516 for ( G4int i = 1 ; i < n1 ; ++i )
517 {
518 if ( E/eV < vE_T[ i ] )
519 {
520 j = i-1;
521 break;
522 }
523 }
524
525 G4double rand_for_mu = G4UniformRand();
526
527 G4int k = 0;
528 for ( G4int i = 0 ; i <= j ; ++i )
529 {
530 G4double Pi = vp_T[ i ] / vp_T[ j ];
531 if ( rand_for_mu < Pi )
532 {
533 k = i;
534 break;
535 }
536 }
537
538 G4double Ei = vE_T[ k ];
539
540 G4double mu = 1 - 2 * Ei / (E/eV) ;
541
542 if ( mu < -1.0 ) mu = -1.0;
543
545 G4double phi = CLHEP::twopi*G4UniformRand();
546 G4double sint= std::sqrt ( 1 - mu*mu );
547 theParticleChange.SetMomentumChange( sint*std::cos(phi), sint*std::sin(phi), mu );
548 }
549 else
550 {
551 // InCoherent Elastic
552
553 // T_L and T_H
554 std::vector<G4double> v_temp;
555 v_temp.clear();
556 for (auto it = incoherentFSs->find(ielement)->second->cbegin();
557 it != incoherentFSs->find(ielement)->second->cend(); ++it)
558 {
559 v_temp.push_back( it->first );
560 }
561
562 // T_L T_H
563 std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
564
565 //
566 // For T_L anEPM_TL and T_H anEPM_TH
567 //
568
569 E_isoAng anEPM_TL_E;
570 E_isoAng anEPM_TH_E;
571
572 if ( tempLH.first != 0.0 && tempLH.second != 0.0 ) {
573 //Interpolate TL and TH
574 anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs->find( ielement )->second->find ( tempLH.first/kelvin )->second );
575 anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs->find( ielement )->second->find ( tempLH.second/kelvin )->second );
576 } else if ( tempLH.first == 0.0 ) {
577 //Extrapolate T0 and T1
578 anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs->find( ielement )->second->find ( v_temp[ 0 ] )->second );
579 anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs->find( ielement )->second->find ( v_temp[ 1 ] )->second );
580 tempLH.first = tempLH.second;
581 tempLH.second = v_temp[ 1 ];
582 } else if ( tempLH.second == 0.0 ) {
583 //Extrapolate Tmax-1 and Tmax
584 anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs->find( ielement )->second->find ( v_temp.back() )->second );
585 auto itv = v_temp.cend();
586 --itv;
587 --itv;
588 anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs->find( ielement )->second->find ( *itv )->second );
589 tempLH.second = tempLH.first;
590 tempLH.first = *itv;
591 }
592
593 // E_isoAng for aTemp and aTrack.GetKineticEnergy()
594 G4double mu=1.0;
595
596 // New Geant4 method - Stochastic interpolation of the final state
597 E_isoAng anEPM_T_E_sampled;
598 G4double rand_temp = G4UniformRand();
599 if ( rand_temp < (aTemp-tempLH.first)/(tempLH.second - tempLH.first) )
600 anEPM_T_E_sampled = anEPM_TH_E;
601 else
602 anEPM_T_E_sampled = anEPM_TL_E;
603
604 mu = getMu ( &anEPM_T_E_sampled );
605
606 // Set Final State
607 theParticleChange.SetEnergyChange( aTrack.GetKineticEnergy() ); // No energy change in Elastic
608 G4double phi = CLHEP::twopi*G4UniformRand();
609 G4double sint= std::sqrt ( 1 - mu*mu );
610 theParticleChange.SetMomentumChange( sint*std::cos(phi), sint*std::sin(phi), mu );
611 }
612 delete dp;
613
614 return &theParticleChange;
615 }
616 else
617 {
618 // Not thermal element
619 // Neutron HP will handle
620 return theHPElastic -> ApplyYourself( aTrack, aNucleus, 1); // L. Thulliez 2021/05/04 (CEA-Saclay)
621 }
622}
623
624
625//**********************************************************
626// Geant4 new algorithm
627//**********************************************************
628
629//--------------------------------------------------
630// New method added by L. Thulliez 2021 (CEA-Saclay)
631//--------------------------------------------------
632std::pair< G4double , G4int> G4ParticleHPThermalScattering::
633sample_inelastic_E( G4double rndm1, G4double rndm2, E_P_E_isoAng* anE_P_E_isoAng )
634{
635 G4int i=0;
636 G4double sE_value=0;
637
638 for ( ; i < anE_P_E_isoAng->secondary_energy_cdf_size-1 ; ++i )
639 {
640 if ( rndm1 >= anE_P_E_isoAng->secondary_energy_cdf[i] &&
641 rndm1 < anE_P_E_isoAng->secondary_energy_cdf[i+1] )
642 {
643 G4double sE_value_i = anE_P_E_isoAng->secondary_energy_value[i];
644 G4double sE_pdf_i = anE_P_E_isoAng->secondary_energy_pdf[i];
645 G4double sE_value_i1 = anE_P_E_isoAng->secondary_energy_value[i+1];
646 G4double sE_pdf_i1 = anE_P_E_isoAng->secondary_energy_pdf[i+1];
647
648 G4double lambda = 0;
649 G4double alpha = (sE_pdf_i1 - sE_pdf_i) / (sE_pdf_i1 + sE_pdf_i);
650 G4double rndm = rndm1;
651
652 if ( std::fabs(alpha) < 1E-8 )
653 {
654 lambda = rndm2;
655 }
656 else
657 {
658 G4double beta = 2 * sE_pdf_i / (sE_pdf_i1 + sE_pdf_i);
659 rndm = rndm2;
660 G4double gamma = -rndm;
661 G4double delta = beta*beta - 4*alpha*gamma;
662
663 if ( delta < 0 && std::fabs(delta) < 1.E-8 ) delta = 0;
664
665 lambda = -beta + std::sqrt(delta);
666 lambda = lambda/(2 * alpha);
667
668 if ( lambda > 1 ) lambda = 1;
669 else if ( lambda < 0 ) lambda = 0;
670 }
671
672 sE_value = sE_value_i + lambda * (sE_value_i1 - sE_value_i);
673
674 break;
675 }
676 }
677
678 return std::pair< G4double , G4int >( sE_value , i );
679}
680
681
682//--------------------------------------------------
683// New method added by L. Thulliez 2021 (CEA-Saclay)
684//--------------------------------------------------
685std::pair< G4double , G4double > G4ParticleHPThermalScattering::
686sample_inelastic_E_mu( G4double pE , std::vector< E_P_E_isoAng* >* vNEP_EPM )
687{
688 // Sample primary energy bin
689 std::map< G4double , G4int > map_energy;
690 map_energy.clear();
691 std::vector< G4double > v_energy;
692 v_energy.clear();
693 G4int i = 0;
694 for (auto itv = vNEP_EPM->cbegin(); itv != vNEP_EPM->cend(); ++itv)
695 {
696 v_energy.push_back( (*itv)->energy );
697 map_energy.insert( std::pair< G4double , G4int >( (*itv)->energy , i ) );
698 i++;
699 }
700
701 std::pair< G4double , G4double > energyLH = find_LH( pE , &v_energy );
702
703 std::vector< E_P_E_isoAng* > pE_P_E_isoAng_limit(2, nullptr);
704
705 if ( energyLH.first != 0.0 && energyLH.second != 0.0 )
706 {
707 pE_P_E_isoAng_limit[0] = (*vNEP_EPM)[ map_energy.find ( energyLH.first )->second ];
708 pE_P_E_isoAng_limit[1] = (*vNEP_EPM)[ map_energy.find ( energyLH.second )->second ];
709 }
710 else if ( energyLH.first == 0.0 )
711 {
712 pE_P_E_isoAng_limit[0] = (*vNEP_EPM)[ 0 ];
713 pE_P_E_isoAng_limit[1] = (*vNEP_EPM)[ 1 ];
714 }
715 if ( energyLH.second == 0.0 )
716 {
717 pE_P_E_isoAng_limit[1] = (*vNEP_EPM).back();
718 auto itv = vNEP_EPM->cend();
719 --itv;
720 --itv;
721 pE_P_E_isoAng_limit[0] = *itv;
722 }
723
724 // Compute interpolation factor of the incident neutron energy
725 G4double factor = (energyLH.second - pE) / (energyLH.second - energyLH.first);
726
727 if ( (energyLH.second - pE) <= 0. && std::fabs(pE/energyLH.second - 1) < 1E-11 ) factor = 0.;
728 if ( (energyLH.first - pE) >= 0. && std::fabs(energyLH.first / pE - 1) < 1E-11 ) factor = 1.;
729
730 G4double rndm1 = G4UniformRand();
731 G4double rndm2 = G4UniformRand();
732
733 // Sample secondary neutron energy
734 std::pair< G4double , G4int > sE_lower = sample_inelastic_E( rndm1, rndm2, pE_P_E_isoAng_limit[0] );
735 std::pair< G4double , G4int > sE_upper = sample_inelastic_E( rndm1, rndm2, pE_P_E_isoAng_limit[1] );
736 G4double sE = factor * sE_lower.first + (1 - factor) * sE_upper.first;
737 sE = sE * eV;
738
739 // Sample cosine knowing the secondary neutron energy
740 rndm1 = G4UniformRand();
741 rndm2 = G4UniformRand();
742 G4double mu_lower = getMu( rndm1, rndm2, pE_P_E_isoAng_limit[0]->vE_isoAngle[sE_lower.second] );
743 G4double mu_upper = getMu( rndm1, rndm2, pE_P_E_isoAng_limit[1]->vE_isoAngle[sE_upper.second] );
744 G4double mu = factor * mu_lower + (1 - factor) * mu_upper;
745
746 return std::pair< G4double , G4double >( sE , mu );
747}
748
749
750//--------------------------------------------------
751// New method added by L. Thulliez 2021 (CEA-Saclay)
752//--------------------------------------------------
753G4double G4ParticleHPThermalScattering::getMu( G4double rndm1, G4double rndm2, E_isoAng* anEPM )
754{
755 G4double result = 0.0;
756
757 G4int in = G4int ( rndm1 * ( (*anEPM).n ) );
758
759 if ( in != 0 )
760 {
761 G4double mu_l = (*anEPM).isoAngle[ in-1 ];
762 G4double mu_h = (*anEPM).isoAngle[ in ];
763 result = ( mu_h - mu_l ) * ( rndm1*((*anEPM).n) - in ) + mu_l;
764 }
765 else
766 {
767 G4double x = rndm1 * (*anEPM).n;
768 G4double ratio = 0.5;
769 if ( x <= ratio )
770 {
771 G4double mu_l = -1;
772 G4double mu_h = (*anEPM).isoAngle[ 0 ];
773 result = ( mu_h - mu_l ) * rndm2 + mu_l;
774 }
775 else
776 {
777 G4double mu_l = (*anEPM).isoAngle[ (*anEPM).n - 1 ];
778 G4double mu_h = 1;
779 result = ( mu_h - mu_l ) * rndm2 + mu_l;
780 }
781 }
782
783 return result;
784}
785
786
787//**********************************************************
788// Geant4 previous algorithm
789//**********************************************************
790
791G4double G4ParticleHPThermalScattering::getMu( E_isoAng* anEPM )
792{
793
794 G4double random = G4UniformRand();
795 G4double result = 0.0;
796
797 G4int in = G4int ( random * ( (*anEPM).n ) );
798
799 if ( in != 0 )
800 {
801 G4double mu_l = (*anEPM).isoAngle[ in-1 ];
802 G4double mu_h = (*anEPM).isoAngle[ in ];
803 result = ( mu_h - mu_l ) * ( random * ( (*anEPM).n ) - in ) + mu_l;
804 }
805 else
806 {
807 G4double x = random * (*anEPM).n;
808 //Bugzilla 1971
809 G4double ratio = 0.5;
810 G4double xx = G4UniformRand();
811 if ( x <= ratio )
812 {
813 G4double mu_l = -1;
814 G4double mu_h = (*anEPM).isoAngle[ 0 ];
815 result = ( mu_h - mu_l ) * xx + mu_l;
816 }
817 else
818 {
819 G4double mu_l = (*anEPM).isoAngle[ (*anEPM).n - 1 ];
820 G4double mu_h = 1;
821 result = ( mu_h - mu_l ) * xx + mu_l;
822 }
823 }
824 return result;
825}
826
827
828std::pair < G4double , G4double > G4ParticleHPThermalScattering::find_LH ( G4double x , std::vector< G4double >* aVector )
829{
830 G4double LL = 0.0;
831 G4double H = 0.0;
832
833 // v->size() == 1 --> LL=H=v(0)
834 if ( aVector->size() == 1 ) {
835 LL = aVector->front();
836 H = aVector->front();
837 } else {
838 // 1) temp < v(0) -> LL=0.0 H=v(0)
839 // 2) v(i-1) < temp <= v(i) -> LL=v(i-1) H=v(i)
840 // 3) v(imax) < temp -> LL=v(imax) H=0.0
841 for ( auto it = aVector->cbegin() ; it != aVector->cend() ; ++it ) {
842 if ( x <= *it ) {
843 H = *it;
844 if ( it != aVector->cbegin() ) {
845 // 2)
846 it--;
847 LL = *it;
848 } else {
849 // 1)
850 LL = 0.0;
851 }
852 break;
853 }
854 }
855 // 3)
856 if ( H == 0.0 ) LL = aVector->back();
857 }
858
859 return std::pair < G4double , G4double > ( LL , H );
860}
861
862
863G4double G4ParticleHPThermalScattering::get_linear_interpolated ( G4double x , std::pair< G4double , G4double > Low , std::pair< G4double , G4double > High )
864{
865 G4double y=0.0;
866 if ( High.first - Low.first != 0 ) {
867 y = ( High.second - Low.second ) / ( High.first - Low.first ) * ( x - Low.first ) + Low.second;
868 } else {
869 if ( High.second == Low.second ) {
870 y = High.second;
871 } else {
872 G4cout << "G4ParticleHPThermalScattering liner interpolation err!!" << G4endl;
873 }
874 }
875
876 return y;
877}
878
879
881G4ParticleHPThermalScattering::create_E_isoAng_from_energy(G4double energy,
882 std::vector<E_isoAng*>* vEPM)
883{
884 E_isoAng anEPM_T_E;
885
886 std::vector<G4double> v_e;
887 v_e.clear();
888 for (auto iv = vEPM->cbegin(); iv != vEPM->cend(); ++iv)
889 v_e.push_back( (*iv)->energy );
890
891 std::pair<G4double, G4double> energyLH = find_LH(energy, &v_e);
892 //G4cout << " " << energy/eV << " " << energyLH.first/eV << " " << energyLH.second/eV << G4endl;
893
894 E_isoAng* panEPM_T_EL = 0;
895 E_isoAng* panEPM_T_EH = 0;
896
897 if (energyLH.first != 0.0 && energyLH.second != 0.0) {
898 for (auto iv = vEPM->cbegin(); iv != vEPM->cend(); ++iv) {
899 if (energyLH.first == (*iv)->energy) {
900 panEPM_T_EL = *iv;
901 ++iv;
902 panEPM_T_EH = *iv;
903 break;
904 }
905 }
906
907 } else if (energyLH.first == 0.0) {
908 panEPM_T_EL = (*vEPM)[0];
909 panEPM_T_EH = (*vEPM)[1];
910
911 } else if (energyLH.second == 0.0) {
912 panEPM_T_EH = (*vEPM).back();
913 auto iv = vEPM->cend();
914 --iv;
915 --iv;
916 panEPM_T_EL = *iv;
917 }
918
919 if (panEPM_T_EL != 0 && panEPM_T_EH != 0) {
920 //checking isoAng has proper values or not
921 // Inelastic/FS, the first and last entries of *vEPM has all zero values.
922 if ( !(check_E_isoAng(panEPM_T_EL) ) ) panEPM_T_EL = panEPM_T_EH;
923 if ( !(check_E_isoAng(panEPM_T_EH) ) ) panEPM_T_EH = panEPM_T_EL;
924
925 if (panEPM_T_EL->n == panEPM_T_EH->n) {
926 anEPM_T_E.energy = energy;
927 anEPM_T_E.n = panEPM_T_EL->n;
928
929 for (G4int i=0; i < panEPM_T_EL->n; ++i) {
930 G4double angle;
931 angle = get_linear_interpolated(energy, std::pair<G4double,G4double>(energyLH.first, panEPM_T_EL->isoAngle[i] ),
932 std::pair<G4double,G4double>(energyLH.second, panEPM_T_EH->isoAngle[i] ) );
933 anEPM_T_E.isoAngle.push_back(angle);
934 }
935
936 } else {
937 G4Exception("G4ParticleHPThermalScattering::create_E_isoAng_from_energy",
938 "NotSupported", JustWarning,
939 "G4ParticleHPThermalScattering does not support yet EL->n != EH->n.");
940 }
941
942 } else {
943 G4Exception("G4ParticleHPThermalScattering::create_E_isoAng_from_energy",
944 "HAD_THERM_000", FatalException,
945 "Pointer panEPM_T_EL or panEPM_T_EH is zero");
946 }
947
948 return anEPM_T_E;
949}
950
951
952G4double G4ParticleHPThermalScattering::
953get_secondary_energy_from_E_P_E_isoAng ( G4double random , E_P_E_isoAng* anE_P_E_isoAng )
954{
955 G4double secondary_energy = 0.0;
956
957 G4int n = anE_P_E_isoAng->n;
958 G4double sum_p = 0.0; // sum_p_H
959 G4double sum_p_L = 0.0;
960
961 G4double total=0.0;
962
963/*
964 delete for speed up
965 for ( G4int i = 0 ; i < n-1 ; ++i )
966 {
967 G4double E_L = anE_P_E_isoAng->vE_isoAngle[i]->energy/eV;
968 G4double E_H = anE_P_E_isoAng->vE_isoAngle[i+1]->energy/eV;
969 G4double dE = E_H - E_L;
970 total += ( ( anE_P_E_isoAng->prob[i] ) * dE );
971 }
972
973 if ( std::abs( total - anE_P_E_isoAng->sum_of_probXdEs ) > 1.0e-14 ) G4cout << total - anE_P_E_isoAng->sum_of_probXdEs << G4endl;
974*/
975 total = anE_P_E_isoAng->sum_of_probXdEs;
976
977 for ( G4int i = 0 ; i < n-1 ; ++i )
978 {
979 G4double E_L = anE_P_E_isoAng->vE_isoAngle[i]->energy/eV;
980 G4double E_H = anE_P_E_isoAng->vE_isoAngle[i+1]->energy/eV;
981 G4double dE = E_H - E_L;
982 sum_p += ( ( anE_P_E_isoAng->prob[i] ) * dE );
983
984 if ( random <= sum_p/total )
985 {
986 secondary_energy = get_linear_interpolated ( random , std::pair < G4double , G4double > ( sum_p_L/total , E_L ) , std::pair < G4double , G4double > ( sum_p/total , E_H ) );
987 secondary_energy = secondary_energy*eV; //need eV
988 break;
989 }
990 sum_p_L = sum_p;
991 }
992
993 return secondary_energy;
994}
995
996
997std::pair< G4double , E_isoAng > G4ParticleHPThermalScattering::
998create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( G4double rand_for_sE , G4double pE , std::vector < E_P_E_isoAng* >* vNEP_EPM )
999{
1000 std::map< G4double , G4int > map_energy;
1001 map_energy.clear();
1002 std::vector< G4double > v_energy;
1003 v_energy.clear();
1004 G4int i = 0;
1005 for (auto itv = vNEP_EPM->cbegin(); itv != vNEP_EPM->cend(); ++itv)
1006 {
1007 v_energy.push_back( (*itv)->energy );
1008 map_energy.insert( std::pair < G4double , G4int > ( (*itv)->energy , i ) );
1009 i++;
1010 }
1011
1012 std::pair < G4double , G4double > energyLH = find_LH ( pE , &v_energy );
1013
1014 E_P_E_isoAng* pE_P_E_isoAng_EL = 0;
1015 E_P_E_isoAng* pE_P_E_isoAng_EH = 0;
1016
1017 if ( energyLH.first != 0.0 && energyLH.second != 0.0 )
1018 {
1019 pE_P_E_isoAng_EL = (*vNEP_EPM)[ map_energy.find ( energyLH.first )->second ];
1020 pE_P_E_isoAng_EH = (*vNEP_EPM)[ map_energy.find ( energyLH.second )->second ];
1021 }
1022 else if ( energyLH.first == 0.0 )
1023 {
1024 pE_P_E_isoAng_EL = (*vNEP_EPM)[ 0 ];
1025 pE_P_E_isoAng_EH = (*vNEP_EPM)[ 1 ];
1026 }
1027 if ( energyLH.second == 0.0 )
1028 {
1029 pE_P_E_isoAng_EH = (*vNEP_EPM).back();
1030 auto itv = vNEP_EPM->cend();
1031 --itv;
1032 --itv;
1033 pE_P_E_isoAng_EL = *itv;
1034 }
1035
1036 G4double sE;
1037 G4double sE_L;
1038 G4double sE_H;
1039
1040 sE_L = get_secondary_energy_from_E_P_E_isoAng ( rand_for_sE , pE_P_E_isoAng_EL );
1041 sE_H = get_secondary_energy_from_E_P_E_isoAng ( rand_for_sE , pE_P_E_isoAng_EH );
1042
1043 sE = get_linear_interpolated ( pE , std::pair < G4double , G4double > ( energyLH.first , sE_L ) , std::pair < G4double , G4double > ( energyLH.second , sE_H ) );
1044
1045
1046 E_isoAng E_isoAng_L = create_E_isoAng_from_energy ( sE , &(pE_P_E_isoAng_EL->vE_isoAngle) );
1047 E_isoAng E_isoAng_H = create_E_isoAng_from_energy ( sE , &(pE_P_E_isoAng_EH->vE_isoAngle) );
1048
1049 E_isoAng anE_isoAng;
1050 //For defeating warning message from compiler
1051 anE_isoAng.n = 1;
1052 anE_isoAng.energy = sE; //never used
1053 if ( E_isoAng_L.n == E_isoAng_H.n )
1054 {
1055 anE_isoAng.n = E_isoAng_L.n;
1056 for ( G4int j=0 ; j < anE_isoAng.n ; ++j )
1057 {
1058 G4double angle;
1059 angle = get_linear_interpolated ( sE , std::pair< G4double , G4double > ( sE_L , E_isoAng_L.isoAngle[ j ] ) , std::pair< G4double , G4double > ( sE_H , E_isoAng_H.isoAngle[ j ] ) );
1060 anE_isoAng.isoAngle.push_back( angle );
1061 }
1062 }
1063 else
1064 {
1065 throw G4HadronicException(__FILE__, __LINE__, "Unexpected values!");
1066 }
1067
1068 return std::pair< G4double , E_isoAng >( sE , anE_isoAng);
1069}
1070
1071
1072void G4ParticleHPThermalScattering::buildPhysicsTable()
1073{
1074 //Is rebuild of physics table a necessity
1075 if ( nMaterial == G4Material::GetMaterialTable()->size() && nElement == G4Element::GetElementTable()->size() ) {
1076 return;
1077 } else {
1078 nMaterial = G4Material::GetMaterialTable()->size();
1079 nElement = G4Element::GetElementTable()->size();
1080 }
1081
1082 dic.clear();
1083 std::map < G4String , G4int > co_dic;
1084
1085 //Searching Nist Materials
1086 static G4ThreadLocal G4MaterialTable* theMaterialTable = nullptr ;
1087 if (!theMaterialTable) theMaterialTable= G4Material::GetMaterialTable();
1088 std::size_t numberOfMaterials = G4Material::GetNumberOfMaterials();
1089 for ( std::size_t i = 0 ; i < numberOfMaterials ; ++i )
1090 {
1091 G4Material* material = (*theMaterialTable)[i];
1092 G4int numberOfElements = (G4int)material->GetNumberOfElements();
1093 for ( G4int j = 0 ; j < numberOfElements ; ++j )
1094 {
1095 const G4Element* element = material->GetElement(j);
1096 if ( names.IsThisThermalElement ( material->GetName() , element->GetName() ) )
1097 {
1098 G4int ts_ID_of_this_geometry;
1099 G4String ts_ndl_name = names.GetTS_NDL_Name( material->GetName() , element->GetName() );
1100 if ( co_dic.find ( ts_ndl_name ) != co_dic.cend() )
1101 {
1102 ts_ID_of_this_geometry = co_dic.find ( ts_ndl_name ) -> second;
1103 }
1104 else
1105 {
1106 ts_ID_of_this_geometry = (G4int)co_dic.size();
1107 co_dic.insert ( std::pair< G4String , G4int >( ts_ndl_name , ts_ID_of_this_geometry ) );
1108 }
1109
1110 //G4cout << "Neutron HP Thermal Scattering: Registering a material-element pair of "
1111 // << material->GetName() << " " << element->GetName()
1112 // << " as internal thermal scattering id of " << ts_ID_of_this_geometry << "." << G4endl;
1113
1114 dic.insert( std::pair < std::pair < G4Material* , const G4Element* > , G4int > ( std::pair < G4Material* , const G4Element* > ( material , element ) , ts_ID_of_this_geometry ) );
1115 }
1116 }
1117 }
1118
1119 //Searching TS Elements
1120 static G4ThreadLocal G4ElementTable* theElementTable = nullptr ;
1121 if (!theElementTable) theElementTable= G4Element::GetElementTable();
1122 std::size_t numberOfElements = G4Element::GetNumberOfElements();
1123 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
1124 {
1125 const G4Element* element = (*theElementTable)[i];
1126 if ( names.IsThisThermalElement ( element->GetName() ) )
1127 {
1128 if ( names.IsThisThermalElement ( element->GetName() ) )
1129 {
1130 G4int ts_ID_of_this_geometry;
1131 G4String ts_ndl_name = names.GetTS_NDL_Name( element->GetName() );
1132 if ( co_dic.find ( ts_ndl_name ) != co_dic.cend() )
1133 {
1134 ts_ID_of_this_geometry = co_dic.find ( ts_ndl_name ) -> second;
1135 }
1136 else
1137 {
1138 ts_ID_of_this_geometry = (G4int)co_dic.size();
1139 co_dic.insert ( std::pair< G4String , G4int >( ts_ndl_name , ts_ID_of_this_geometry ) );
1140 }
1141 dic.insert( std::pair < std::pair < const G4Material* , const G4Element* > , G4int > ( std::pair < const G4Material* , const G4Element* > ( (G4Material*)nullptr , element ) , ts_ID_of_this_geometry ) );
1142 }
1143 }
1144 }
1145
1146 G4cout << G4endl;
1147 G4cout << "Neutron HP Thermal Scattering: Following material-element pairs or elements are registered." << G4endl;
1148 for ( std::map < std::pair < const G4Material* , const G4Element* > , G4int >::iterator it = dic.begin() ; it != dic.end() ; it++ )
1149 {
1150 if ( it->first.first != nullptr )
1151 {
1152 G4cout << "Material " << it->first.first->GetName() << " - Element " << it->first.second->GetName() << ", internal thermal scattering id " << it->second << G4endl;
1153 }
1154 else
1155 {
1156 G4cout << "Element " << it->first.second->GetName() << ", internal thermal scattering id " << it->second << G4endl;
1157 }
1158 }
1159 G4cout << G4endl;
1160
1161 // Read Cross Section Data files
1162
1164 coherentFSs = hpmanager->GetThermalScatteringCoherentFinalStates();
1165 incoherentFSs = hpmanager->GetThermalScatteringIncoherentFinalStates();
1166 inelasticFSs = hpmanager->GetThermalScatteringInelasticFinalStates();
1167
1169
1170 clearCurrentFSData();
1171
1172 if ( coherentFSs == nullptr ) coherentFSs = new std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* >;
1173 if ( incoherentFSs == nullptr ) incoherentFSs = new std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* >;
1174 if ( inelasticFSs == nullptr ) inelasticFSs = new std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* >;
1175
1176 G4String dirName;
1177 if ( !G4FindDataDir( "G4NEUTRONHPDATA" ) )
1178 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
1179 dirName = G4FindDataDir( "G4NEUTRONHPDATA" );
1180
1181 for (auto it = co_dic.cbegin() ; it != co_dic.cend() ; ++it)
1182 {
1183 G4String tsndlName = it->first;
1184 G4int ts_ID = it->second;
1185
1186 // Coherent
1187 G4String fsName = "/ThermalScattering/Coherent/FS/";
1188 G4String fileName = dirName + fsName + tsndlName;
1189 coherentFSs->insert ( std::pair < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* > ( ts_ID , readACoherentFSDATA( fileName ) ) );
1190
1191 // incoherent elastic
1192 fsName = "/ThermalScattering/Incoherent/FS/";
1193 fileName = dirName + fsName + tsndlName;
1194 incoherentFSs->insert ( std::pair < G4int , std::map < G4double , std::vector < E_isoAng* >* >* > ( ts_ID , readAnIncoherentFSDATA( fileName ) ) );
1195
1196 // inelastic
1197 fsName = "/ThermalScattering/Inelastic/FS/";
1198 fileName = dirName + fsName + tsndlName;
1199 inelasticFSs->insert ( std::pair < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* > ( ts_ID , readAnInelasticFSDATA( fileName ) ) );
1200 }
1201
1202 hpmanager->RegisterThermalScatteringCoherentFinalStates( coherentFSs );
1203 hpmanager->RegisterThermalScatteringIncoherentFinalStates( incoherentFSs );
1204 hpmanager->RegisterThermalScatteringInelasticFinalStates( inelasticFSs );
1205 }
1206
1207 theXSection->BuildPhysicsTable( *(G4Neutron::Neutron()) );
1208}
1209
1210
1211G4int G4ParticleHPThermalScattering::getTS_ID ( const G4Material* material , const G4Element* element )
1212{
1213 G4int result = -1;
1214 if ( dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) ) != dic.end() )
1215 result = dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) )->second;
1216 return result;
1217}
1218
1219
1220const std::pair<G4double, G4double> G4ParticleHPThermalScattering::GetFatalEnergyCheckLevels() const
1221{
1222 // max energy non-conservation is mass of heavy nucleus
1223 return std::pair<G4double, G4double>(10.0*perCent, 350.0*CLHEP::GeV);
1224}
1225
1226
1228{
1229 names.AddThermalElement( nameG4Element , filename );
1230 theXSection->AddUserThermalScatteringFile( nameG4Element , filename );
1231 buildPhysicsTable();
1232}
1233
1234
1235G4bool G4ParticleHPThermalScattering::check_E_isoAng( E_isoAng* anE_IsoAng )
1236{
1237 G4bool result=false;
1238
1239 G4int n = anE_IsoAng->n;
1240 G4double sum=0.0;
1241 for ( G4int i = 0 ; i < n ; ++i ) {
1242 sum += anE_IsoAng->isoAngle[ i ];
1243 }
1244 if ( sum != 0.0 ) result = true;
1245
1246 return result;
1247}
1248
1249
1251{
1252 outFile << "High Precision model based on thermal scattering data in\n"
1253 << "evaluated nuclear data libraries for neutrons below 5eV\n"
1254 << "on specific materials\n";
1255}
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
const G4String & GetName() const
Definition: G4Element.hh:127
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
G4double GetTemperature() const
Definition: G4Material.hh:177
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
void BuildPhysicsTable(const G4ParticleDefinition &)
void RegisterThermalScatteringCoherentFinalStates(std::map< G4int, std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * > *val)
void RegisterThermalScatteringIncoherentFinalStates(std::map< G4int, std::map< G4double, std::vector< E_isoAng * > * > * > *val)
void RegisterThermalScatteringInelasticFinalStates(std::map< G4int, std::map< G4double, std::vector< E_P_E_isoAng * > * > * > *val)
std::map< G4int, std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * > * GetThermalScatteringCoherentFinalStates()
static G4ParticleHPManager * GetInstance()
std::map< G4int, std::map< G4double, std::vector< E_isoAng * > * > * > * GetThermalScatteringIncoherentFinalStates()
void GetDataStream(G4String, std::istringstream &iss)
std::map< G4int, std::map< G4double, std::vector< E_P_E_isoAng * > * > * > * GetThermalScatteringInelasticFinalStates()
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetCoherentCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
G4double GetInelasticCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
void AddUserThermalScatteringFile(G4String, G4String)
virtual void ModelDescription(std::ostream &outFile) const
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double total(Particle const *const p1, Particle const *const p2)
G4double energy(const ThreeVector &p, const G4double m)
G4bool IsMasterThread()
Definition: G4Threading.cc:124
std::vector< G4double > prob
std::vector< E_isoAng * > vE_isoAngle
std::vector< G4double > secondary_energy_pdf
std::vector< G4double > secondary_energy_cdf
std::vector< G4double > secondary_energy_value
std::vector< G4double > isoAngle
#define G4ThreadLocal
Definition: tls.hh:77