Geant4  10.01
G4NeutronHPInelastic.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 // this code implementation is the intellectual property of
27 // neutron_hp -- source file
28 // J.P. Wellisch, Nov-1996
29 // A prototype of the low energy neutron transport model.
30 //
31 // By copying, distributing or modifying the Program (or any work
32 // based on the Program) you indicate your acceptance of this statement,
33 // and all its terms.
34 //
35 //
36 // 070523 bug fix for G4FPE_DEBUG on by A. Howard (and T. Koi)
37 // 081203 limit maximum trial for creating final states add protection for 1H isotope case by T. Koi
38 //
39 #include "G4NeutronHPInelastic.hh"
40 #include "G4SystemOfUnits.hh"
41 
42 #include "G4NeutronHPManager.hh"
43 
45  :G4HadronicInteraction("NeutronHPInelastic")
46  ,theInelastic(NULL)
47  ,numEle(0)
48  {
49  SetMinEnergy( 0.0 );
50  SetMaxEnergy( 20.*MeV );
51 /*
52 
53  G4int istatus;
54 #if defined WIN32-VC
55  istatus = system("echo %G4NEUTRONHPDATA%");
56 #else
57  istatus = system("echo $G4NEUTRONHPDATA");
58 #endif
59  if ( istatus < 0 )
60  {
61  G4cout << "Warning! system(\"echo $G4NEUTRONHPDATA\") returns error value at G4NeutronHPInelastic" << G4endl;
62  }
63 
64 // G4cout << " entering G4NeutronHPInelastic constructor"<<G4endl;
65  if(!getenv("G4NEUTRONHPDATA"))
66  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
67  dirName = getenv("G4NEUTRONHPDATA");
68  G4String tString = "/Inelastic";
69  dirName = dirName + tString;
70  numEle = G4Element::GetNumberOfElements();
71 */
72 /*
73  theInelastic = new G4NeutronHPChannelList[numEle];
74  for (G4int i=0; i<numEle; i++)
75  {
76  theInelastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
77  G4int itry = 0;
78  do
79  {
80  theInelastic[i].Register(&theNFS, "F01"); // has
81  theInelastic[i].Register(&theNXFS, "F02");
82  theInelastic[i].Register(&the2NDFS, "F03");
83  theInelastic[i].Register(&the2NFS, "F04"); // has, E Done
84  theInelastic[i].Register(&the3NFS, "F05"); // has, E Done
85  theInelastic[i].Register(&theNAFS, "F06");
86  theInelastic[i].Register(&theN3AFS, "F07");
87  theInelastic[i].Register(&the2NAFS, "F08");
88  theInelastic[i].Register(&the3NAFS, "F09");
89  theInelastic[i].Register(&theNPFS, "F10");
90  theInelastic[i].Register(&theN2AFS, "F11");
91  theInelastic[i].Register(&the2N2AFS, "F12");
92  theInelastic[i].Register(&theNDFS, "F13");
93  theInelastic[i].Register(&theNTFS, "F14");
94  theInelastic[i].Register(&theNHe3FS, "F15");
95  theInelastic[i].Register(&theND2AFS, "F16");
96  theInelastic[i].Register(&theNT2AFS, "F17");
97  theInelastic[i].Register(&the4NFS, "F18"); // has, E Done
98  theInelastic[i].Register(&the2NPFS, "F19");
99  theInelastic[i].Register(&the3NPFS, "F20");
100  theInelastic[i].Register(&theN2PFS, "F21");
101  theInelastic[i].Register(&theNPAFS, "F22");
102  theInelastic[i].Register(&thePFS, "F23");
103  theInelastic[i].Register(&theDFS, "F24");
104  theInelastic[i].Register(&theTFS, "F25");
105  theInelastic[i].Register(&theHe3FS, "F26");
106  theInelastic[i].Register(&theAFS, "F27");
107  theInelastic[i].Register(&the2AFS, "F28");
108  theInelastic[i].Register(&the3AFS, "F29");
109  theInelastic[i].Register(&the2PFS, "F30");
110  theInelastic[i].Register(&thePAFS, "F31");
111  theInelastic[i].Register(&theD2AFS, "F32");
112  theInelastic[i].Register(&theT2AFS, "F33");
113  theInelastic[i].Register(&thePDFS, "F34");
114  theInelastic[i].Register(&thePTFS, "F35");
115  theInelastic[i].Register(&theDAFS, "F36");
116  theInelastic[i].RestartRegistration();
117  itry++;
118  }
119  //while(!theInelastic[i].HasDataInAnyFinalState());
120  while( !theInelastic[i].HasDataInAnyFinalState() && itry < 6 );
121  // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK
122 
123  if ( itry == 6 )
124  {
125  // No Final State at all.
126  G4bool exceptional = false;
127  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
128  {
129  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
130  }
131  if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
132  }
133  }
134 */
135 /*
136 
137  for (G4int i=0; i<numEle; i++)
138  {
139  theInelastic.push_back( new G4NeutronHPChannelList );
140  (*theInelastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
141  G4int itry = 0;
142  do
143  {
144  (*theInelastic[i]).Register(&theNFS, "F01"); // has
145  (*theInelastic[i]).Register(&theNXFS, "F02");
146  (*theInelastic[i]).Register(&the2NDFS, "F03");
147  (*theInelastic[i]).Register(&the2NFS, "F04"); // has, E Done
148  (*theInelastic[i]).Register(&the3NFS, "F05"); // has, E Done
149  (*theInelastic[i]).Register(&theNAFS, "F06");
150  (*theInelastic[i]).Register(&theN3AFS, "F07");
151  (*theInelastic[i]).Register(&the2NAFS, "F08");
152  (*theInelastic[i]).Register(&the3NAFS, "F09");
153  (*theInelastic[i]).Register(&theNPFS, "F10");
154  (*theInelastic[i]).Register(&theN2AFS, "F11");
155  (*theInelastic[i]).Register(&the2N2AFS, "F12");
156  (*theInelastic[i]).Register(&theNDFS, "F13");
157  (*theInelastic[i]).Register(&theNTFS, "F14");
158  (*theInelastic[i]).Register(&theNHe3FS, "F15");
159  (*theInelastic[i]).Register(&theND2AFS, "F16");
160  (*theInelastic[i]).Register(&theNT2AFS, "F17");
161  (*theInelastic[i]).Register(&the4NFS, "F18"); // has, E Done
162  (*theInelastic[i]).Register(&the2NPFS, "F19");
163  (*theInelastic[i]).Register(&the3NPFS, "F20");
164  (*theInelastic[i]).Register(&theN2PFS, "F21");
165  (*theInelastic[i]).Register(&theNPAFS, "F22");
166  (*theInelastic[i]).Register(&thePFS, "F23");
167  (*theInelastic[i]).Register(&theDFS, "F24");
168  (*theInelastic[i]).Register(&theTFS, "F25");
169  (*theInelastic[i]).Register(&theHe3FS, "F26");
170  (*theInelastic[i]).Register(&theAFS, "F27");
171  (*theInelastic[i]).Register(&the2AFS, "F28");
172  (*theInelastic[i]).Register(&the3AFS, "F29");
173  (*theInelastic[i]).Register(&the2PFS, "F30");
174  (*theInelastic[i]).Register(&thePAFS, "F31");
175  (*theInelastic[i]).Register(&theD2AFS, "F32");
176  (*theInelastic[i]).Register(&theT2AFS, "F33");
177  (*theInelastic[i]).Register(&thePDFS, "F34");
178  (*theInelastic[i]).Register(&thePTFS, "F35");
179  (*theInelastic[i]).Register(&theDAFS, "F36");
180  (*theInelastic[i]).RestartRegistration();
181  itry++;
182  }
183  while( !(*theInelastic[i]).HasDataInAnyFinalState() && itry < 6 );
184  // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK
185 
186  if ( itry == 6 )
187  {
188  // No Final State at all.
189  G4bool exceptional = false;
190  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
191  {
192  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
193  }
194  if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
195  }
196 
197  }
198 */
199  }
200 
202  {
203 // delete [] theInelastic;
204  for ( std::vector<G4NeutronHPChannelList*>::iterator
205  it = (*theInelastic).begin() ; it != (*theInelastic).end() ; it++ )
206  {
207  delete *it;
208  }
209  (*theInelastic).clear();
210  }
211 
212  #include "G4NeutronHPThermalBoost.hh"
213 
215  {
216  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
218  const G4Material * theMaterial = aTrack.GetMaterial();
219  G4int n = theMaterial->GetNumberOfElements();
220  G4int index = theMaterial->GetElement(0)->GetIndex();
221  G4int it=0;
222  if(n!=1)
223  {
224  xSec = new G4double[n];
225  G4double sum=0;
226  G4int i;
227  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
228  G4double rWeight;
229  G4NeutronHPThermalBoost aThermalE;
230  for (i=0; i<n; i++)
231  {
232  index = theMaterial->GetElement(i)->GetIndex();
233  rWeight = NumAtomsPerVolume[i];
234  //xSec[i] = theInelastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
235  xSec[i] = (*(*theInelastic)[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
236  theMaterial->GetElement(i),
237  theMaterial->GetTemperature()));
238  xSec[i] *= rWeight;
239  sum+=xSec[i];
240  }
241  G4double random = G4UniformRand();
242  G4double running = 0;
243  for (i=0; i<n; i++)
244  {
245  running += xSec[i];
246  index = theMaterial->GetElement(i)->GetIndex();
247  it = i;
248  //if(random<=running/sum) break;
249  if( sum == 0 || random<=running/sum) break;
250  }
251  delete [] xSec;
252  }
253 
254  //return theInelastic[index].ApplyYourself(theMaterial->GetElement(it), aTrack);
255  G4HadFinalState* result = (*(*theInelastic)[index]).ApplyYourself(theMaterial->GetElement(it), aTrack);
256 
257  //Overwrite target parameters
258  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
259  const G4Element* target_element = (*G4Element::GetElementTable())[index];
260  const G4Isotope* target_isotope=NULL;
261  G4int iele = target_element->GetNumberOfIsotopes();
262  for ( G4int j = 0 ; j != iele ; j++ ) {
263  target_isotope=target_element->GetIsotope( j );
264  if ( target_isotope->GetN() == G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
265  }
266  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
267  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
268  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
269  aNucleus.SetIsotope( target_isotope );
270 
272  return result;
273  }
274 
275 const std::pair<G4double, G4double> G4NeutronHPInelastic::GetFatalEnergyCheckLevels() const
276 {
277  // max energy non-conservation is mass of heavy nucleus
278 // if ( getenv("G4NEUTRONHP_DO_NOT_ADJUST_FINAL_STATE") ) return std::pair<G4double, G4double>(5*perCent,250*GeV);
279  // This should be same to the hadron default value
280 // return std::pair<G4double, G4double>(10*perCent,10*GeV);
281  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
282 }
283 
284 /*
285 void G4NeutronHPInelastic::addChannelForNewElement()
286 {
287  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
288  {
289  G4cout << "G4NeutronHPInelastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
290 
291  theInelastic.push_back( new G4NeutronHPChannelList );
292  (*theInelastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
293  G4int itry = 0;
294  do
295  {
296  (*theInelastic[i]).Register(&theNFS, "F01"); // has
297  (*theInelastic[i]).Register(&theNXFS, "F02");
298  (*theInelastic[i]).Register(&the2NDFS, "F03");
299  (*theInelastic[i]).Register(&the2NFS, "F04"); // has, E Done
300  (*theInelastic[i]).Register(&the3NFS, "F05"); // has, E Done
301  (*theInelastic[i]).Register(&theNAFS, "F06");
302  (*theInelastic[i]).Register(&theN3AFS, "F07");
303  (*theInelastic[i]).Register(&the2NAFS, "F08");
304  (*theInelastic[i]).Register(&the3NAFS, "F09");
305  (*theInelastic[i]).Register(&theNPFS, "F10");
306  (*theInelastic[i]).Register(&theN2AFS, "F11");
307  (*theInelastic[i]).Register(&the2N2AFS, "F12");
308  (*theInelastic[i]).Register(&theNDFS, "F13");
309  (*theInelastic[i]).Register(&theNTFS, "F14");
310  (*theInelastic[i]).Register(&theNHe3FS, "F15");
311  (*theInelastic[i]).Register(&theND2AFS, "F16");
312  (*theInelastic[i]).Register(&theNT2AFS, "F17");
313  (*theInelastic[i]).Register(&the4NFS, "F18"); // has, E Done
314  (*theInelastic[i]).Register(&the2NPFS, "F19");
315  (*theInelastic[i]).Register(&the3NPFS, "F20");
316  (*theInelastic[i]).Register(&theN2PFS, "F21");
317  (*theInelastic[i]).Register(&theNPAFS, "F22");
318  (*theInelastic[i]).Register(&thePFS, "F23");
319  (*theInelastic[i]).Register(&theDFS, "F24");
320  (*theInelastic[i]).Register(&theTFS, "F25");
321  (*theInelastic[i]).Register(&theHe3FS, "F26");
322  (*theInelastic[i]).Register(&theAFS, "F27");
323  (*theInelastic[i]).Register(&the2AFS, "F28");
324  (*theInelastic[i]).Register(&the3AFS, "F29");
325  (*theInelastic[i]).Register(&the2PFS, "F30");
326  (*theInelastic[i]).Register(&thePAFS, "F31");
327  (*theInelastic[i]).Register(&theD2AFS, "F32");
328  (*theInelastic[i]).Register(&theT2AFS, "F33");
329  (*theInelastic[i]).Register(&thePDFS, "F34");
330  (*theInelastic[i]).Register(&thePTFS, "F35");
331  (*theInelastic[i]).Register(&theDAFS, "F36");
332  (*theInelastic[i]).RestartRegistration();
333  itry++;
334  }
335  while( !(*theInelastic[i]).HasDataInAnyFinalState() && itry < 6 );
336  // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK
337 
338  if ( itry == 6 )
339  {
340  // No Final State at all.
341  G4bool exceptional = false;
342  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
343  {
344  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
345  }
346  if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
347  }
348  }
349 
350  numEle = (G4int)G4Element::GetNumberOfElements();
351 }
352 */
353 
355 {
357 }
359 {
361 }
363 {
365 
366  theInelastic = hpmanager->GetInelasticFinalStates();
367 
368  if ( !G4Threading::IsWorkerThread() ) {
369 
370  if ( theInelastic == NULL ) theInelastic = new std::vector<G4NeutronHPChannelList*>;
371 
372  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
373 
374  if ( theInelastic->size() == G4Element::GetNumberOfElements() ) {
376  return;
377  }
378 
379  if ( !getenv("G4NEUTRONHPDATA") )
380  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
381  dirName = getenv("G4NEUTRONHPDATA");
382  G4String tString = "/Inelastic";
383  dirName = dirName + tString;
384 
385  for (G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++)
386  {
387  (*theInelastic).push_back( new G4NeutronHPChannelList );
388  (*(*theInelastic)[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
389  G4int itry = 0;
390  do
391  {
392  (*(*theInelastic)[i]).Register(&theNFS, "F01"); // has
393  (*(*theInelastic)[i]).Register(&theNXFS, "F02");
394  (*(*theInelastic)[i]).Register(&the2NDFS, "F03");
395  (*(*theInelastic)[i]).Register(&the2NFS, "F04"); // has, E Done
396  (*(*theInelastic)[i]).Register(&the3NFS, "F05"); // has, E Done
397  (*(*theInelastic)[i]).Register(&theNAFS, "F06");
398  (*(*theInelastic)[i]).Register(&theN3AFS, "F07");
399  (*(*theInelastic)[i]).Register(&the2NAFS, "F08");
400  (*(*theInelastic)[i]).Register(&the3NAFS, "F09");
401  (*(*theInelastic)[i]).Register(&theNPFS, "F10");
402  (*(*theInelastic)[i]).Register(&theN2AFS, "F11");
403  (*(*theInelastic)[i]).Register(&the2N2AFS, "F12");
404  (*(*theInelastic)[i]).Register(&theNDFS, "F13");
405  (*(*theInelastic)[i]).Register(&theNTFS, "F14");
406  (*(*theInelastic)[i]).Register(&theNHe3FS, "F15");
407  (*(*theInelastic)[i]).Register(&theND2AFS, "F16");
408  (*(*theInelastic)[i]).Register(&theNT2AFS, "F17");
409  (*(*theInelastic)[i]).Register(&the4NFS, "F18"); // has, E Done
410  (*(*theInelastic)[i]).Register(&the2NPFS, "F19");
411  (*(*theInelastic)[i]).Register(&the3NPFS, "F20");
412  (*(*theInelastic)[i]).Register(&theN2PFS, "F21");
413  (*(*theInelastic)[i]).Register(&theNPAFS, "F22");
414  (*(*theInelastic)[i]).Register(&thePFS, "F23");
415  (*(*theInelastic)[i]).Register(&theDFS, "F24");
416  (*(*theInelastic)[i]).Register(&theTFS, "F25");
417  (*(*theInelastic)[i]).Register(&theHe3FS, "F26");
418  (*(*theInelastic)[i]).Register(&theAFS, "F27");
419  (*(*theInelastic)[i]).Register(&the2AFS, "F28");
420  (*(*theInelastic)[i]).Register(&the3AFS, "F29");
421  (*(*theInelastic)[i]).Register(&the2PFS, "F30");
422  (*(*theInelastic)[i]).Register(&thePAFS, "F31");
423  (*(*theInelastic)[i]).Register(&theD2AFS, "F32");
424  (*(*theInelastic)[i]).Register(&theT2AFS, "F33");
425  (*(*theInelastic)[i]).Register(&thePDFS, "F34");
426  (*(*theInelastic)[i]).Register(&thePTFS, "F35");
427  (*(*theInelastic)[i]).Register(&theDAFS, "F36");
428  (*(*theInelastic)[i]).RestartRegistration();
429  itry++;
430  }
431  while( !(*(*theInelastic)[i]).HasDataInAnyFinalState() && itry < 6 );
432  // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK
433 
434  if ( itry == 6 )
435  {
436  // No Final State at all.
437  G4bool exceptional = false;
438  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
439  {
440  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
441  }
442  if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
443  }
444 
445  }
447 
448  }
450 }
G4NeutronHPNDInelasticFS theNDFS
G4NeutronHPN2PInelasticFS theN2PFS
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4NeutronHP2PInelasticFS the2PFS
void RegisterInelasticFinalStates(std::vector< G4NeutronHPChannelList * > *val)
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
static const double MeV
Definition: G4SIunits.hh:193
G4NeutronHPHe3InelasticFS theHe3FS
G4NeutronHPPTInelasticFS thePTFS
G4NeutronHPT2AInelasticFS theT2AFS
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4NeutronHP2N2AInelasticFS the2N2AFS
G4NeutronHPDAInelasticFS theDAFS
G4NeutronHPNPAInelasticFS theNPAFS
std::vector< G4NeutronHPChannelList * > * theInelastic
G4NeutronHPNInelasticFS theNFS
G4NeutronHPTInelasticFS theTFS
G4NeutronHPPDInelasticFS thePDFS
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4NeutronHP2NAInelasticFS the2NAFS
G4NeutronHPNXInelasticFS theNXFS
G4NeutronHP2NPInelasticFS the2NPFS
G4NeutronHPD2AInelasticFS theD2AFS
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
void BuildPhysicsTable(const G4ParticleDefinition &)
G4NeutronHP3NInelasticFS the3NFS
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4NeutronHP4NInelasticFS the4NFS
G4NeutronHP2NDInelasticFS the2NDFS
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4NeutronHP3NAInelasticFS the3NAFS
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:95
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
G4NeutronHPNT2AInelasticFS theNT2AFS
G4NeutronHPPInelasticFS thePFS
bool G4bool
Definition: G4Types.hh:79
G4NeutronHPND2AInelasticFS theND2AFS
G4NeutronHP2AInelasticFS the2AFS
G4NeutronHPNAInelasticFS theNAFS
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
G4NeutronHPDInelasticFS theDFS
G4bool IsWorkerThread()
Definition: G4Threading.cc:128
void SetVerboseLevel(G4int i)
G4NeutronHPNHe3InelasticFS theNHe3FS
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4NeutronHPPAInelasticFS thePAFS
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4NeutronHP3AInelasticFS the3AFS
G4double GetTemperature() const
Definition: G4Material.hh:180
G4NeutronHPN2AInelasticFS theN2AFS
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
G4NeutronHPNPInelasticFS theNPFS
G4NeutronHP3NPInelasticFS the3NPFS
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4NeutronHPN3AInelasticFS theN3AFS
#define DBL_MAX
Definition: templates.hh:83
G4NeutronHPNTInelasticFS theNTFS
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
G4NeutronHP2NInelasticFS the2NFS
G4NeutronHPAInelasticFS theAFS
std::vector< G4NeutronHPChannelList * > * GetInelasticFinalStates()