Geant4  10.01.p02
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 //
46 #include "G4ParticleHPManager.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "G4Neutron.hh"
49 #include "G4ElementTable.hh"
50 #include "G4MaterialTable.hh"
51 
53  :G4HadronicInteraction("NeutronHPThermalScattering")
54 {
55 
57 
58  SetMinEnergy( 0.*eV );
59  SetMaxEnergy( 4*eV );
62 
65 
66 }
67 
68 
69 
71 {
72 
73  for ( std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* >::iterator it = incoherentFSs.begin() ; it != incoherentFSs.end() ; it++ )
74  {
75  std::map < G4double , std::vector < E_isoAng* >* >::iterator itt;
76  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
77  {
78  std::vector< E_isoAng* >::iterator ittt;
79  for ( ittt = itt->second->begin(); ittt != itt->second->end() ; ittt++ )
80  {
81  delete *ittt;
82  }
83  delete itt->second;
84  }
85  delete it->second;
86  }
87 
88  for ( std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* >::iterator it = coherentFSs.begin() ; it != coherentFSs.end() ; it++ )
89  {
90  std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >::iterator itt;
91  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
92  {
93  std::vector < std::pair< G4double , G4double >* >::iterator ittt;
94  for ( ittt = itt->second->begin(); ittt != itt->second->end() ; ittt++ )
95  {
96  delete *ittt;
97  }
98  delete itt->second;
99  }
100  delete it->second;
101  }
102 
103  for ( std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* >::iterator it = inelasticFSs.begin() ; it != inelasticFSs.end() ; it++ )
104  {
105  std::map < G4double , std::vector < E_P_E_isoAng* >* >::iterator itt;
106  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
107  {
108  std::vector < E_P_E_isoAng* >::iterator ittt;
109  for ( ittt = itt->second->begin(); ittt != itt->second->end() ; ittt++ )
110  {
111  std::vector < E_isoAng* >::iterator it4;
112  for ( it4 = (*ittt)->vE_isoAngle.begin() ; it4 != (*ittt)->vE_isoAngle.end() ; it4++ )
113  {
114  delete *it4;
115  }
116  delete *ittt;
117  }
118  delete itt->second;
119  }
120  delete it->second;
121  }
122 
123  delete theHPElastic;
124  delete theXSection;
125 }
126 
127 
128 
129 std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* G4ParticleHPThermalScattering::readACoherentFSDATA( G4String name )
130 {
131 
132  std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* aCoherentFSDATA = new std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >;
133 
134  //std::ifstream theChannel( name.c_str() );
135  std::istringstream theChannel(std::ios::in);
137 
138  std::vector< G4double > vBraggE;
139 
140  G4int dummy;
141  while ( theChannel >> dummy ) // MF
142  {
143  theChannel >> dummy; // MT
144  G4double temp;
145  theChannel >> temp;
146  std::vector < std::pair< G4double , G4double >* >* anBragE_P = new std::vector < std::pair< G4double , G4double >* >;
147 
148  G4int n;
149  theChannel >> n;
150  for ( G4int i = 0 ; i < n ; i++ )
151  {
152  G4double Ei;
153  G4double Pi;
154  if ( aCoherentFSDATA->size() == 0 )
155  {
156  theChannel >> Ei;
157  vBraggE.push_back( Ei );
158  }
159  else
160  {
161  Ei = vBraggE[ i ];
162  }
163  theChannel >> Pi;
164  anBragE_P->push_back ( new std::pair < G4double , G4double > ( Ei , Pi ) );
165  //G4cout << "Coherent Elastic " << Ei << " " << Pi << G4endl;
166  }
167  aCoherentFSDATA->insert ( std::pair < G4double , std::vector < std::pair< G4double , G4double >* >* > ( temp , anBragE_P ) );
168  }
169 
170  return aCoherentFSDATA;
171 }
172 
173 
174 
175 std::map < G4double , std::vector < E_P_E_isoAng* >* >* G4ParticleHPThermalScattering::readAnInelasticFSDATA ( G4String name )
176 {
177  std::map < G4double , std::vector < E_P_E_isoAng* >* >* anT_E_P_E_isoAng = new std::map < G4double , std::vector < E_P_E_isoAng* >* >;
178 
179  //std::ifstream theChannel( name.c_str() );
180  std::istringstream theChannel(std::ios::in);
182 
183  G4int dummy;
184  while ( theChannel >> dummy ) // MF
185  {
186  theChannel >> dummy; // MT
187  G4double temp;
188  theChannel >> temp;
189  std::vector < E_P_E_isoAng* >* vE_P_E_isoAng = new std::vector < E_P_E_isoAng* >;
190  G4int n;
191  theChannel >> n;
192  for ( G4int i = 0 ; i < n ; i++ )
193  {
194  vE_P_E_isoAng->push_back ( readAnE_P_E_isoAng ( &theChannel ) );
195  }
196  anT_E_P_E_isoAng->insert ( std::pair < G4double , std::vector < E_P_E_isoAng* >* > ( temp , vE_P_E_isoAng ) );
197  }
198  //theChannel.close();
199 
200  return anT_E_P_E_isoAng;
201 }
202 
203 
204 
206 {
207  E_P_E_isoAng* aData = new E_P_E_isoAng;
208 
209  G4double dummy;
211  G4int nep , nl;
212  *file >> dummy;
213  *file >> energy;
214  aData->energy = energy*eV;
215  *file >> dummy;
216  *file >> dummy;
217  *file >> nep;
218  *file >> nl;
219  aData->n = nep/nl;
220  for ( G4int i = 0 ; i < aData->n ; i++ )
221  {
222  G4double prob;
223  E_isoAng* anE_isoAng = new E_isoAng;
224  aData->vE_isoAngle.push_back( anE_isoAng );
225  *file >> energy;
226  anE_isoAng->energy = energy*eV;
227  anE_isoAng->n = nl - 2;
228  anE_isoAng->isoAngle.resize( anE_isoAng->n );
229  *file >> prob;
230  aData->prob.push_back( prob );
231  //G4cout << "G4ParticleHPThermalScattering inelastic " << energy/eV << " " << i << " " << prob << " " << aData->prob[ i ] << G4endl;
232  for ( G4int j = 0 ; j < anE_isoAng->n ; j++ )
233  {
234  G4double x;
235  *file >> x;
236  anE_isoAng->isoAngle[j] = x ;
237  //G4cout << "G4ParticleHPThermalScattering inelastic " << x << anE_isoAng->isoAngle[j] << G4endl;
238  }
239  }
240 
241  // Calcuate sum_of_provXdEs
242  G4double total = 0;
243  for ( G4int i = 0 ; i < aData->n - 1 ; i++ )
244  {
245  G4double E_L = aData->vE_isoAngle[i]->energy/eV;
246  G4double E_H = aData->vE_isoAngle[i+1]->energy/eV;
247  G4double dE = E_H - E_L;
248  total += ( ( aData->prob[i] ) * dE );
249  }
250  aData->sum_of_probXdEs = total;
251 
252  return aData;
253 }
254 
255 
256 
257 std::map < G4double , std::vector < E_isoAng* >* >* G4ParticleHPThermalScattering::readAnIncoherentFSDATA ( G4String name )
258 {
259  std::map < G4double , std::vector < E_isoAng* >* >* T_E = new std::map < G4double , std::vector < E_isoAng* >* >;
260 
261  //std::ifstream theChannel( name.c_str() );
262  std::istringstream theChannel(std::ios::in);
264 
265  G4int dummy;
266  while ( theChannel >> dummy ) // MF
267  {
268  theChannel >> dummy; // MT
269  G4double temp;
270  theChannel >> temp;
271  std::vector < E_isoAng* >* vE_isoAng = new std::vector < E_isoAng* >;
272  G4int n;
273  theChannel >> n;
274  for ( G4int i = 0 ; i < n ; i++ )
275  vE_isoAng->push_back ( readAnE_isoAng( &theChannel ) );
276  T_E->insert ( std::pair < G4double , std::vector < E_isoAng* >* > ( temp , vE_isoAng ) );
277  }
278  //theChannel.close();
279 
280  return T_E;
281 }
282 
283 
284 
286 {
287  E_isoAng* aData = new E_isoAng;
288 
289  G4double dummy;
291  G4int n;
292  *file >> dummy;
293  *file >> energy;
294  *file >> dummy;
295  *file >> dummy;
296  *file >> n;
297  *file >> dummy;
298  aData->energy = energy*eV;
299  aData->n = n-2;
300  aData->isoAngle.resize( n );
301 
302  *file >> dummy;
303  *file >> dummy;
304  for ( G4int i = 0 ; i < aData->n ; i++ )
305  *file >> aData->isoAngle[i];
306 
307  return aData;
308 }
309 
310 
311 
313 {
314 
315  //Trick for dynamically generated materials
320  }
321 // Select Element > Reaction >
322 
323  const G4Material * theMaterial = aTrack.GetMaterial();
324  G4double aTemp = theMaterial->GetTemperature();
325  G4int n = theMaterial->GetNumberOfElements();
326  //static const G4ElementTable* theElementTable = G4Element::GetElementTable();
327 
328  G4bool findThermalElement = false;
329  G4int ielement;
330  const G4Element* theElement = NULL;
331  for ( G4int i = 0; i < n ; i++ )
332  {
333  theElement = theMaterial->GetElement(i);
334  //Select target element
335  if ( aNucleus.GetZ_asInt() == (G4int)(theElement->GetZ() + 0.5 ) )
336  {
337  //Check Applicability of Thermal Scattering
338  if ( getTS_ID( NULL , theElement ) != -1 )
339  {
340  ielement = getTS_ID( NULL , theElement );
341  findThermalElement = true;
342  break;
343  }
344  else if ( getTS_ID( theMaterial , theElement ) != -1 )
345  {
346  ielement = getTS_ID( theMaterial , theElement );
347  findThermalElement = true;
348  break;
349  }
350  }
351  }
352 
353  if ( findThermalElement == true )
354  {
355 
356 // Select Reaction (Inelastic, coherent, incoherent)
357 
358  const G4ParticleDefinition* pd = aTrack.GetDefinition();
359  G4DynamicParticle* dp = new G4DynamicParticle ( pd , aTrack.Get4Momentum() );
360  G4double total = theXSection->GetCrossSection( dp , theElement , theMaterial );
361  G4double inelastic = theXSection->GetInelasticCrossSection( dp , theElement , theMaterial );
362 
363 
364  G4double random = G4UniformRand();
365  if ( random <= inelastic/total )
366  {
367  // Inelastic
368 
369  // T_L and T_H
370  std::map < G4double , std::vector< E_P_E_isoAng* >* >::iterator it;
371  std::vector<G4double> v_temp;
372  v_temp.clear();
373  for ( it = inelasticFSs.find( ielement )->second->begin() ; it != inelasticFSs.find( ielement )->second->end() ; it++ )
374  {
375  v_temp.push_back( it->first );
376  }
377 
378 // T_L T_H
379  std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
380 //
381 // For T_L aNEP_EPM_TL and T_H aNEP_EPM_TH
382 //
383  std::vector< E_P_E_isoAng* >* vNEP_EPM_TL = 0;
384  std::vector< E_P_E_isoAng* >* vNEP_EPM_TH = 0;
385 
386  if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
387  {
388  vNEP_EPM_TL = inelasticFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second;
389  vNEP_EPM_TH = inelasticFSs.find( ielement )->second->find ( tempLH.second/kelvin )->second;
390  }
391  else if ( tempLH.first == 0.0 )
392  {
393  std::map < G4double , std::vector< E_P_E_isoAng* >* >::iterator itm;
394  itm = inelasticFSs.find( ielement )->second->begin();
395  vNEP_EPM_TL = itm->second;
396  itm++;
397  vNEP_EPM_TH = itm->second;
398  tempLH.first = tempLH.second;
399  tempLH.second = itm->first;
400  }
401  else if ( tempLH.second == 0.0 )
402  {
403  std::map < G4double , std::vector< E_P_E_isoAng* >* >::iterator itm;
404  itm = inelasticFSs.find( ielement )->second->end();
405  itm--;
406  vNEP_EPM_TH = itm->second;
407  itm--;
408  vNEP_EPM_TL = itm->second;
409  tempLH.second = tempLH.first;
410  tempLH.first = itm->first;
411  }
412 
413  G4double rand_for_sE = G4UniformRand();
414 
415  std::pair< G4double , E_isoAng > TL = create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( rand_for_sE , aTrack.GetKineticEnergy() , vNEP_EPM_TL );
416  std::pair< G4double , E_isoAng > TH = create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( rand_for_sE , aTrack.GetKineticEnergy() , vNEP_EPM_TH );
417 
418  G4double sE;
419  sE = get_linear_interpolated ( aTemp , std::pair < G4double , G4double > ( tempLH.first , TL.first ) , std::pair < G4double , G4double > ( tempLH.second , TH.first ) );
420 
421  G4double mu=1.0;
422  E_isoAng anE_isoAng;
423  if ( TL.second.n == TH.second.n )
424  {
425  anE_isoAng.energy = sE;
426  anE_isoAng.n = TL.second.n;
427  for ( G4int i=0 ; i < anE_isoAng.n ; i++ )
428  {
429  G4double angle;
430  angle = get_linear_interpolated ( aTemp , std::pair< G4double , G4double > ( tempLH.first , TL.second.isoAngle[ i ] ) , std::pair< G4double , G4double > ( tempLH.second , TH.second.isoAngle[ i ] ) );
431  anE_isoAng.isoAngle.push_back( angle );
432  }
433  mu = getMu( &anE_isoAng );
434 
435  } else {
436  //TL.second.n != TH.second.n
437  G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Do not yet supported");
438  }
439 
440  //set
442  theParticleChange.SetMomentumChange( 0.0 , std::sqrt ( 1 - mu*mu ) , mu );
443 
444  }
445  //else if ( random <= ( inelastic + theXSection->GetCoherentCrossSection( dp , (*theElementTable)[ ielement ] , aTemp ) ) / total )
446  else if ( random <= ( inelastic + theXSection->GetCoherentCrossSection( dp , theElement , theMaterial ) ) / total )
447  {
448  // Coherent Elastic
449 
450  G4double E = aTrack.GetKineticEnergy();
451 
452  // T_L and T_H
453  std::map < G4double , std::vector< std::pair< G4double , G4double >* >* >::iterator it;
454  std::vector<G4double> v_temp;
455  v_temp.clear();
456  for ( it = coherentFSs.find( ielement )->second->begin() ; it != coherentFSs.find( ielement )->second->end() ; it++ )
457  {
458  v_temp.push_back( it->first );
459  }
460 
461 // T_L T_H
462  std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
463 //
464 //
465 // For T_L anEPM_TL and T_H anEPM_TH
466 //
467  std::vector< std::pair< G4double , G4double >* >* pvE_p_TL = NULL;
468  std::vector< std::pair< G4double , G4double >* >* pvE_p_TH = NULL;
469 
470  if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
471  {
472  pvE_p_TL = coherentFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second;
473  pvE_p_TH = coherentFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second;
474  }
475  else if ( tempLH.first == 0.0 )
476  {
477  pvE_p_TL = coherentFSs.find( ielement )->second->find ( v_temp[ 0 ] )->second;
478  pvE_p_TH = coherentFSs.find( ielement )->second->find ( v_temp[ 1 ] )->second;
479  tempLH.first = tempLH.second;
480  tempLH.second = v_temp[ 1 ];
481  }
482  else if ( tempLH.second == 0.0 )
483  {
484  pvE_p_TL = coherentFSs.find( ielement )->second->find ( v_temp.back() )->second;
485  std::vector< G4double >::iterator itv;
486  itv = v_temp.end();
487  itv--;
488  itv--;
489  pvE_p_TL = coherentFSs.find( ielement )->second->find ( *itv )->second;
490  tempLH.second = tempLH.first;
491  tempLH.first = *itv;
492  }
493  else
494  {
495  //tempLH.first == 0.0 && tempLH.second
496  G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Unexpected temperature values in data");
497  }
498 
499  std::vector< G4double > vE_T;
500  std::vector< G4double > vp_T;
501 
502  G4int n1 = pvE_p_TL->size();
503  //G4int n2 = pvE_p_TH->size();
504 
505  for ( G4int i=1 ; i < n1 ; i++ )
506  {
507  if ( (*pvE_p_TL)[i]->first != (*pvE_p_TH)[i]->first ) G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data!");
508  vE_T.push_back ( (*pvE_p_TL)[i]->first );
509  vp_T.push_back ( get_linear_interpolated ( aTemp , std::pair< G4double , G4double > ( tempLH.first , (*pvE_p_TL)[i]->second ) , std::pair< G4double , G4double > ( tempLH.second , (*pvE_p_TL)[i]->second ) ) );
510  }
511 
512  G4int j = 0;
513  for ( G4int i = 1 ; i < n ; i++ )
514  {
515  if ( E/eV < vE_T[ i ] )
516  {
517  j = i-1;
518  break;
519  }
520  }
521 
522  G4double rand_for_mu = G4UniformRand();
523 
524  G4int k = 0;
525  for ( G4int i = 1 ; i < j ; i++ )
526  {
527  G4double Pi = vp_T[ i ] / vp_T[ j ];
528  if ( rand_for_mu < Pi )
529  {
530  k = i-1;
531  break;
532  }
533  }
534 
535  //G4double Ei = vE_T[ j ];
536  G4double Ei = vE_T[ k ];
537 
538  G4double mu = 1 - 2 * Ei / (E/eV) ;
539  //111102
540  if ( mu < -1.0 ) mu = -1.0;
541 
543  theParticleChange.SetMomentumChange( 0.0 , std::sqrt ( 1 - mu*mu ) , mu );
544 
545 
546  }
547  else
548  {
549  // InCoherent Elastic
550 
551  // T_L and T_H
552  std::map < G4double , std::vector < E_isoAng* >* >::iterator it;
553  std::vector<G4double> v_temp;
554  v_temp.clear();
555  for ( it = incoherentFSs.find( ielement )->second->begin() ; it != incoherentFSs.find( ielement )->second->end() ; it++ )
556  {
557  v_temp.push_back( it->first );
558  }
559 
560 // T_L T_H
561  std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
562 
563 //
564 // For T_L anEPM_TL and T_H anEPM_TH
565 //
566 
567  E_isoAng anEPM_TL_E;
568  E_isoAng anEPM_TH_E;
569 
570  if ( tempLH.first != 0.0 && tempLH.second != 0.0 ) {
571  //Interpolate TL and TH
572  anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second );
573  anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( tempLH.second/kelvin )->second );
574  } else if ( tempLH.first == 0.0 ) {
575  //Extrapolate T0 and T1
576  anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( v_temp[ 0 ] )->second );
577  anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( v_temp[ 1 ] )->second );
578  tempLH.first = tempLH.second;
579  tempLH.second = v_temp[ 1 ];
580  } else if ( tempLH.second == 0.0 ) {
581  //Extrapolate Tmax-1 and Tmax
582  anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( v_temp.back() )->second );
583  std::vector< G4double >::iterator itv;
584  itv = v_temp.end();
585  itv--;
586  itv--;
587  anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( *itv )->second );
588  tempLH.second = tempLH.first;
589  tempLH.first = *itv;
590  }
591 
592  // E_isoAng for aTemp and aTrack.GetKineticEnergy()
593  G4double mu=1.0;
594  E_isoAng anEPM_T_E;
595 
596  if ( anEPM_TL_E.n == anEPM_TH_E.n )
597  {
598  anEPM_T_E.n = anEPM_TL_E.n;
599  for ( G4int i=0 ; i < anEPM_TL_E.n ; i++ )
600  {
601  G4double angle;
602  angle = get_linear_interpolated ( aTemp , std::pair< G4double , G4double > ( tempLH.first , anEPM_TL_E.isoAngle[ i ] ) , std::pair< G4double , G4double > ( tempLH.second , anEPM_TH_E.isoAngle[ i ] ) );
603  anEPM_T_E.isoAngle.push_back( angle );
604  }
605  mu = getMu ( &anEPM_T_E );
606 
607  } else {
608  // anEPM_TL_E.n != anEPM_TH_E.n
609  G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Do not yet supported");
610  }
611 
612  // Set Final State
613  theParticleChange.SetEnergyChange( aTrack.GetKineticEnergy() ); // No energy change in Elastic
614  theParticleChange.SetMomentumChange( 0.0 , std::sqrt ( 1 - mu*mu ) , mu );
615 
616  }
617  delete dp;
618 
619  return &theParticleChange;
620 
621  }
622  else
623  {
624  // Not thermal element
625  // Neutron HP will handle
626  return theHPElastic -> ApplyYourself( aTrack, aNucleus );
627  }
628 
629 }
630 
631 
632 
634 {
635 
636  G4double random = G4UniformRand();
637  G4double result = 0.0;
638 
639  G4int in = int ( random * ( (*anEPM).n ) );
640 
641  if ( in != 0 )
642  {
643  G4double mu_l = (*anEPM).isoAngle[ in-1 ];
644  G4double mu_h = (*anEPM).isoAngle[ in ];
645  result = ( mu_h - mu_l ) * ( random * ( (*anEPM).n ) - in ) + mu_l;
646  }
647  else
648  {
649  G4double x = random * (*anEPM).n;
650  G4double D = ( (*anEPM).isoAngle[ 0 ] - ( -1 ) ) + ( 1 - (*anEPM).isoAngle[ (*anEPM).n - 1 ] );
651  G4double ratio = ( (*anEPM).isoAngle[ 0 ] - ( -1 ) ) / D;
652  if ( x <= ratio )
653  {
654  G4double mu_l = -1;
655  G4double mu_h = (*anEPM).isoAngle[ 0 ];
656  result = ( mu_h - mu_l ) * x + mu_l;
657  }
658  else
659  {
660  G4double mu_l = (*anEPM).isoAngle[ (*anEPM).n - 1 ];
661  G4double mu_h = 1;
662  result = ( mu_h - mu_l ) * x + mu_l;
663  }
664  }
665  return result;
666 }
667 
668 
669 
670 std::pair < G4double , G4double > G4ParticleHPThermalScattering::find_LH ( G4double x , std::vector< G4double >* aVector )
671 {
672  G4double L = 0.0;
673  G4double H = 0.0;
674 
675  // v->size() == 1 --> L=H=v(0)
676  if ( aVector->size() == 1 ) {
677  L = aVector->front();
678  H = aVector->front();
679  } else {
680  // 1) temp < v(0) -> L=0.0 H=v(0)
681  // 2) v(i-1) < temp <= v(i) -> L=v(i-1) H=v(i)
682  // 3) v(imax) < temp -> L=v(imax) H=0.0
683  for ( std::vector< G4double >::iterator
684  it = aVector->begin() ; it != aVector->end() ; it++ ) {
685  if ( x <= *it ) {
686  H = *it;
687  if ( it != aVector->begin() ) {
688  // 2)
689  it--;
690  L = *it;
691  } else {
692  // 1)
693  L = 0.0;
694  }
695  break;
696  }
697  }
698  // 3)
699  if ( H == 0.0 ) L = aVector->back();
700  }
701 
702  return std::pair < G4double , G4double > ( L , H );
703 }
704 
705 
706 
707 G4double G4ParticleHPThermalScattering::get_linear_interpolated ( G4double x , std::pair< G4double , G4double > Low , std::pair< G4double , G4double > High )
708 {
709  G4double y=0.0;
710  if ( High.first - Low.first != 0 ) {
711  y = ( High.second - Low.second ) / ( High.first - Low.first ) * ( x - Low.first ) + Low.second;
712  } else {
713  if ( High.second == Low.second ) {
714  y = High.second;
715  } else {
716  G4cout << "G4ParticleHPThermalScattering liner interpolation err!!" << G4endl;
717  }
718  }
719 
720  return y;
721 }
722 
723 
724 
726 {
727  E_isoAng anEPM_T_E;
728 
729  std::vector< E_isoAng* >::iterator iv;
730 
731  std::vector< G4double > v_e;
732  v_e.clear();
733  for ( iv = vEPM->begin() ; iv != vEPM->end() ; iv++ )
734  v_e.push_back ( (*iv)->energy );
735 
736  std::pair < G4double , G4double > energyLH = find_LH ( energy , &v_e );
737  //G4cout << " " << energy/eV << " " << energyLH.first/eV << " " << energyLH.second/eV << G4endl;
738 
739  E_isoAng* panEPM_T_EL=0;
740  E_isoAng* panEPM_T_EH=0;
741 
742  if ( energyLH.first != 0.0 && energyLH.second != 0.0 )
743  {
744  for ( iv = vEPM->begin() ; iv != vEPM->end() ; iv++ )
745  {
746  if ( energyLH.first == (*iv)->energy )
747  break;
748  }
749  panEPM_T_EL = *iv;
750  iv++;
751  panEPM_T_EH = *iv;
752  }
753  else if ( energyLH.first == 0.0 )
754  {
755  panEPM_T_EL = (*vEPM)[0];
756  panEPM_T_EH = (*vEPM)[1];
757  }
758  else if ( energyLH.second == 0.0 )
759  {
760  panEPM_T_EH = (*vEPM).back();
761  iv = vEPM->end();
762  iv--;
763  iv--;
764  panEPM_T_EL = *iv;
765  }
766 
767  if ( panEPM_T_EL->n == panEPM_T_EH->n )
768  {
769  anEPM_T_E.energy = energy;
770  anEPM_T_E.n = panEPM_T_EL->n;
771 
772  for ( G4int i=0 ; i < panEPM_T_EL->n ; i++ )
773  {
774  G4double angle;
775  angle = get_linear_interpolated ( energy , std::pair< G4double , G4double > ( energyLH.first , panEPM_T_EL->isoAngle[ i ] ) , std::pair< G4double , G4double > ( energyLH.second , panEPM_T_EH->isoAngle[ i ] ) );
776  anEPM_T_E.isoAngle.push_back( angle );
777  }
778  }
779  else
780  {
781  G4cout << "G4ParticleHPThermalScattering Do not Suuport yet." << G4endl;
782  }
783 
784 
785  return anEPM_T_E;
786 }
787 
788 
789 
791 {
792 
793  G4double secondary_energy = 0.0;
794 
795  G4int n = anE_P_E_isoAng->n;
796  G4double sum_p = 0.0; // sum_p_H
797  G4double sum_p_L = 0.0;
798 
799  G4double total=0.0;
800 
801 /*
802  delete for speed up
803  for ( G4int i = 0 ; i < n-1 ; i++ )
804  {
805  G4double E_L = anE_P_E_isoAng->vE_isoAngle[i]->energy/eV;
806  G4double E_H = anE_P_E_isoAng->vE_isoAngle[i+1]->energy/eV;
807  G4double dE = E_H - E_L;
808  total += ( ( anE_P_E_isoAng->prob[i] ) * dE );
809  }
810 
811  if ( std::abs( total - anE_P_E_isoAng->sum_of_probXdEs ) > 1.0e-14 ) G4cout << total - anE_P_E_isoAng->sum_of_probXdEs << G4endl;
812 */
813  total = anE_P_E_isoAng->sum_of_probXdEs;
814 
815  for ( G4int i = 0 ; i < n-1 ; i++ )
816  {
817  G4double E_L = anE_P_E_isoAng->vE_isoAngle[i]->energy/eV;
818  G4double E_H = anE_P_E_isoAng->vE_isoAngle[i+1]->energy/eV;
819  G4double dE = E_H - E_L;
820  sum_p += ( ( anE_P_E_isoAng->prob[i] ) * dE );
821 
822  if ( random <= sum_p/total )
823  {
824  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 ) );
825  secondary_energy = secondary_energy*eV; //need eV
826  break;
827  }
828  sum_p_L = sum_p;
829  }
830 
831  return secondary_energy;
832 }
833 
834 
835 
836 std::pair< G4double , E_isoAng > G4ParticleHPThermalScattering::create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( G4double rand_for_sE , G4double pE , std::vector < E_P_E_isoAng* >* vNEP_EPM )
837 {
838 
839  std::map< G4double , G4int > map_energy;
840  map_energy.clear();
841  std::vector< G4double > v_energy;
842  v_energy.clear();
843  std::vector< E_P_E_isoAng* >::iterator itv;
844  G4int i = 0;
845  for ( itv = vNEP_EPM->begin(); itv != vNEP_EPM->end(); itv++ )
846  {
847  v_energy.push_back( (*itv)->energy );
848  map_energy.insert( std::pair < G4double , G4int > ( (*itv)->energy , i ) );
849  i++;
850  }
851 
852  std::pair < G4double , G4double > energyLH = find_LH ( pE , &v_energy );
853 
854  E_P_E_isoAng* pE_P_E_isoAng_EL = 0;
855  E_P_E_isoAng* pE_P_E_isoAng_EH = 0;
856 
857  if ( energyLH.first != 0.0 && energyLH.second != 0.0 )
858  {
859  pE_P_E_isoAng_EL = (*vNEP_EPM)[ map_energy.find ( energyLH.first )->second ];
860  pE_P_E_isoAng_EH = (*vNEP_EPM)[ map_energy.find ( energyLH.second )->second ];
861  }
862  else if ( energyLH.first == 0.0 )
863  {
864  pE_P_E_isoAng_EL = (*vNEP_EPM)[ 0 ];
865  pE_P_E_isoAng_EH = (*vNEP_EPM)[ 1 ];
866  }
867  if ( energyLH.second == 0.0 )
868  {
869  pE_P_E_isoAng_EH = (*vNEP_EPM).back();
870  itv = vNEP_EPM->end();
871  itv--;
872  itv--;
873  pE_P_E_isoAng_EL = *itv;
874  }
875 
876 
877  G4double sE;
878  G4double sE_L;
879  G4double sE_H;
880 
881 
882  sE_L = get_secondary_energy_from_E_P_E_isoAng ( rand_for_sE , pE_P_E_isoAng_EL );
883  sE_H = get_secondary_energy_from_E_P_E_isoAng ( rand_for_sE , pE_P_E_isoAng_EH );
884 
885  sE = get_linear_interpolated ( pE , std::pair < G4double , G4double > ( energyLH.first , sE_L ) , std::pair < G4double , G4double > ( energyLH.second , sE_H ) );
886 
887 
888  E_isoAng E_isoAng_L = create_E_isoAng_from_energy ( sE , &(pE_P_E_isoAng_EL->vE_isoAngle) );
889  E_isoAng E_isoAng_H = create_E_isoAng_from_energy ( sE , &(pE_P_E_isoAng_EH->vE_isoAngle) );
890 
891  E_isoAng anE_isoAng;
892  //For defeating warning message from compiler
893  anE_isoAng.n = 1;
894  anE_isoAng.energy = sE; //never used
895  if ( E_isoAng_L.n == E_isoAng_H.n )
896  {
897  anE_isoAng.n = E_isoAng_L.n;
898  for ( G4int j=0 ; j < anE_isoAng.n ; j++ )
899  {
900  G4double angle;
901  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 ] ) );
902  anE_isoAng.isoAngle.push_back( angle );
903  }
904  }
905  else
906  {
907  //G4cout << "Do not Suuport yet." << G4endl;
908  throw G4HadronicException(__FILE__, __LINE__, "Unexpected values!");
909  }
910 
911 
912 
913  return std::pair< G4double , E_isoAng >( sE , anE_isoAng);
914 }
915 
917 {
918 
919  dic.clear();
920  std::map < G4String , G4int > co_dic;
921 
922  //Searching Nist Materials
923  static G4ThreadLocal G4MaterialTable* theMaterialTable = 0 ; if (!theMaterialTable) theMaterialTable= G4Material::GetMaterialTable();
924  size_t numberOfMaterials = G4Material::GetNumberOfMaterials();
925  for ( size_t i = 0 ; i < numberOfMaterials ; i++ )
926  {
927  G4Material* material = (*theMaterialTable)[i];
928  size_t numberOfElements = material->GetNumberOfElements();
929  for ( size_t j = 0 ; j < numberOfElements ; j++ )
930  {
931  const G4Element* element = material->GetElement(j);
932  if ( names.IsThisThermalElement ( material->GetName() , element->GetName() ) )
933  {
934  G4int ts_ID_of_this_geometry;
935  G4String ts_ndl_name = names.GetTS_NDL_Name( material->GetName() , element->GetName() );
936  if ( co_dic.find ( ts_ndl_name ) != co_dic.end() )
937  {
938  ts_ID_of_this_geometry = co_dic.find ( ts_ndl_name ) -> second;
939  }
940  else
941  {
942  ts_ID_of_this_geometry = co_dic.size();
943  co_dic.insert ( std::pair< G4String , G4int >( ts_ndl_name , ts_ID_of_this_geometry ) );
944  }
945 
946  //G4cout << "Neutron HP Thermal Scattering: Registering a material-element pair of "
947  // << material->GetName() << " " << element->GetName()
948  // << " as internal thermal scattering id of " << ts_ID_of_this_geometry << "." << G4endl;
949 
950  dic.insert( std::pair < std::pair < G4Material* , const G4Element* > , G4int > ( std::pair < G4Material* , const G4Element* > ( material , element ) , ts_ID_of_this_geometry ) );
951  }
952  }
953  }
954 
955  //Searching TS Elements
956  static G4ThreadLocal G4ElementTable* theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
957  size_t numberOfElements = G4Element::GetNumberOfElements();
958  //size_t numberOfThermalElements = 0;
959  for ( size_t i = 0 ; i < numberOfElements ; i++ )
960  {
961  const G4Element* element = (*theElementTable)[i];
962  if ( names.IsThisThermalElement ( element->GetName() ) )
963  {
964  if ( names.IsThisThermalElement ( element->GetName() ) )
965  {
966  G4int ts_ID_of_this_geometry;
967  G4String ts_ndl_name = names.GetTS_NDL_Name( element->GetName() );
968  if ( co_dic.find ( ts_ndl_name ) != co_dic.end() )
969  {
970  ts_ID_of_this_geometry = co_dic.find ( ts_ndl_name ) -> second;
971  }
972  else
973  {
974  ts_ID_of_this_geometry = co_dic.size();
975  co_dic.insert ( std::pair< G4String , G4int >( ts_ndl_name , ts_ID_of_this_geometry ) );
976  }
977 
978  //G4cout << "Neutron HP Thermal Scattering: Registering an element of "
979  // << material->GetName() << " " << element->GetName()
980  // << " as internal thermal scattering id of " << ts_ID_of_this_geometry << "." << G4endl;
981 
982  dic.insert( std::pair < std::pair < const G4Material* , const G4Element* > , G4int > ( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , element ) , ts_ID_of_this_geometry ) );
983  }
984  }
985  }
986 
987  G4cout << G4endl;
988  G4cout << "Neutron HP Thermal Scattering: Following material-element pairs or elements are registered." << G4endl;
989  for ( std::map < std::pair < const G4Material* , const G4Element* > , G4int >::iterator it = dic.begin() ; it != dic.end() ; it++ )
990  {
991  if ( it->first.first != NULL )
992  {
993  G4cout << "Material " << it->first.first->GetName() << " - Element " << it->first.second->GetName() << ", internal thermal scattering id " << it->second << G4endl;
994  }
995  else
996  {
997  G4cout << "Element " << it->first.second->GetName() << ", internal thermal scattering id " << it->second << G4endl;
998  }
999  }
1000  G4cout << G4endl;
1001 
1002  // Read Cross Section Data files
1003 
1004  G4String dirName;
1005  if ( !getenv( "G4NEUTRONHPDATA" ) )
1006  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
1007  dirName = getenv( "G4NEUTRONHPDATA" );
1008 
1009  //dirName = baseName + "/ThermalScattering";
1010 
1011  G4String name;
1012 
1013  for ( std::map < G4String , G4int >::iterator it = co_dic.begin() ; it != co_dic.end() ; it++ )
1014  {
1015  G4String tsndlName = it->first;
1016  G4int ts_ID = it->second;
1017 
1018  // Coherent
1019  G4String fsName = "/ThermalScattering/Coherent/FS/";
1020  G4String fileName = dirName + fsName + tsndlName;
1021  coherentFSs.insert ( std::pair < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* > ( ts_ID , readACoherentFSDATA( fileName ) ) );
1022 
1023  // incoherent elastic
1024  fsName = "/ThermalScattering/Incoherent/FS/";
1025  fileName = dirName + fsName + tsndlName;
1026  incoherentFSs.insert ( std::pair < G4int , std::map < G4double , std::vector < E_isoAng* >* >* > ( ts_ID , readAnIncoherentFSDATA( fileName ) ) );
1027 
1028  // inelastic
1029  fsName = "/ThermalScattering/Inelastic/FS/";
1030  fileName = dirName + fsName + tsndlName;
1031  inelasticFSs.insert ( std::pair < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* > ( ts_ID , readAnInelasticFSDATA( fileName ) ) );
1032  }
1033 
1035 }
1036 
1037 
1039 {
1040  G4int result = -1;
1041  if ( dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) ) != dic.end() )
1042  result = dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) )->second;
1043  return result;
1044 }
1045 
1046 const std::pair<G4double, G4double> G4ParticleHPThermalScattering::GetFatalEnergyCheckLevels() const
1047 {
1048  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
1049  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
1050 }
1051 
1053 {
1054  names.AddThermalElement( nameG4Element , filename );
1055  theXSection->AddUserThermalScatteringFile( nameG4Element , filename );
1057 }
static G4ParticleHPManager * GetInstance()
G4ParticleHPThermalScatteringNames names
std::vector< G4double > isoAngle
std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * readACoherentFSDATA(G4String)
G4int first(char) const
G4double get_secondary_energy_from_E_P_E_isoAng(G4double, E_P_E_isoAng *)
std::map< G4int, std::map< G4double, std::vector< E_isoAng * > * > * > incoherentFSs
G4String name
Definition: TRTMaterials.hh:40
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetZ() const
Definition: G4Element.hh:131
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:588
std::vector< G4Material * > G4MaterialTable
void GetDataStream(G4String, std::istringstream &iss)
E_isoAng create_E_isoAng_from_energy(G4double, std::vector< E_isoAng * > *)
#define G4ThreadLocal
Definition: tls.hh:89
std::pair< G4double, E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng(G4double, G4double, std::vector< E_P_E_isoAng * > *)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
std::map< G4double, std::vector< E_P_E_isoAng * > * > * readAnInelasticFSDATA(G4String)
G4double GetCoherentCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
static const G4int L[nN]
void SetMinEnergy(G4double anEnergy)
static const G4double dE
std::map< G4double, std::vector< E_isoAng * > * > * readAnIncoherentFSDATA(G4String)
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
const G4ParticleDefinition * GetDefinition() const
std::map< G4int, std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * > coherentFSs
E_P_E_isoAng * readAnE_P_E_isoAng(std::istream *)
G4ParticleHPThermalScatteringData * theXSection
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
G4int getTS_ID(const G4Material *, const G4Element *)
static const double second
Definition: G4SIunits.hh:138
static const double perCent
Definition: G4SIunits.hh:296
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4int n
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:595
const G4LorentzVector & Get4Momentum() const
static const double kelvin
Definition: G4SIunits.hh:260
static const double eV
Definition: G4SIunits.hh:194
void SetEnergyChange(G4double anEnergy)
G4double total(Particle const *const p1, Particle const *const p2)
G4double get_linear_interpolated(G4double, std::pair< G4double, G4double >, std::pair< G4double, G4double >)
void AddUserThermalScatteringFile(G4String, G4String)
G4double energy(const ThreeVector &p, const G4double m)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
std::map< G4int, std::map< G4double, std::vector< E_P_E_isoAng * > * > * > inelasticFSs
std::pair< G4double, G4double > find_LH(G4double, std::vector< G4double > *)
void SetMaxEnergy(const G4double anEnergy)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
G4double GetInelasticCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
const G4String & GetName() const
Definition: G4Element.hh:127
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
#define DBL_MAX
Definition: templates.hh:83
std::map< std::pair< const G4Material *, const G4Element * >, G4int > dic
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
void BuildPhysicsTable(const G4ParticleDefinition &)
void SetMomentumChange(const G4ThreeVector &aV)
std::vector< E_isoAng * > vE_isoAngle
std::vector< G4double > prob