Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPInelasticBaseFS.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // neutron_hp -- source file
27 // J.P. Wellisch, Nov-1996
28 // A prototype of the low energy neutron transport model.
29 //
30 // 080801 Give a warning message for irregular mass value in data file by T. Koi
31 // Introduce theNDLDataA,Z which has A and Z of NDL data by T. Koi
32 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
33 // 101111 Add Special treatment for Be9(n,2n)Be8(2a) case by T. Koi
34 //
36 #include "G4NeutronHPManager.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4Nucleus.hh"
39 #include "G4NucleiProperties.hh"
40 #include "G4He3.hh"
41 #include "G4Alpha.hh"
42 #include "G4Electron.hh"
43 #include "G4NeutronHPDataUsed.hh"
44 
45 #include "G4IonTable.hh"
46 
48 {
49  // char the[100] = {""};
50  // std::ostrstream ost(the, 100, std::ios::out);
51  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
52  // G4String * aName = new G4String(the);
53  // std::ifstream from(*aName, std::ios::in);
54 
55  std::ostringstream ost;
56  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
57  G4String aName = ost.str();
58  std::ifstream from(aName, std::ios::in);
59 
60  if(!from) return; // no data found for this isotope
61  // std::ifstream theGammaData(*aName, std::ios::in);
62  std::ifstream theGammaData(aName, std::ios::in);
63 
64  G4double eps = 0.001;
66  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) -
67  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
68  theGammas.Init(theGammaData);
69  // delete aName;
70 }
71 
73 {
74  gammaPath = "/Inelastic/Gammas/";
75  if(!getenv("G4NEUTRONHPDATA"))
76  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
77  G4String tBase = getenv("G4NEUTRONHPDATA");
78  gammaPath = tBase+gammaPath;
79  G4String tString = dirName;
80  G4bool dbool;
81  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M,tString, bit, dbool);
82  G4String filename = aFile.GetName();
83  SetAZMs( A, Z, M, aFile);
84  //theBaseA = aFile.GetA();
85  //theBaseZ = aFile.GetZ();
86  // theNDLDataA = (int)aFile.GetA();
87  // theNDLDataZ = aFile.GetZ();
88  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
89  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
90  {
91  if(getenv("NeutronHPNamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
92  hasAnyData = false;
93  hasFSData = false;
94  hasXsec = false;
95  return;
96  }
97  //theBaseA = A;
98  //theBaseZ = G4int(Z+.5);
99  //std::ifstream theData(filename, std::ios::in);
100  //130205 For compressed data files
101  std::istringstream theData(std::ios::in);
102  G4NeutronHPManager::GetInstance()->GetDataStream(filename,theData);
103  //130205 END
104  if(!theData) //"!" is a operator of ios
105  {
106  hasAnyData = false;
107  hasFSData = false;
108  hasXsec = false;
109  //theData.close();
110  return; // no data for exactly this isotope and FS
111  }
112  // here we go
113  G4int infoType, dataType, dummy=INT_MAX;
114  hasFSData = false;
115  while (theData >> infoType)
116  {
117  theData >> dataType;
118  if(dummy==INT_MAX) theData >> dummy >> dummy;
119  if(dataType==3)
120  {
121  G4int total;
122  theData >> total;
123  theXsection->Init(theData, total, eV);
124  }
125  else if(dataType==4)
126  {
128  theAngularDistribution->Init(theData);
129  hasFSData = true;
130  }
131  else if(dataType==5)
132  {
134  theEnergyDistribution->Init(theData);
135  hasFSData = true;
136  }
137  else if(dataType==6)
138  {
140  theEnergyAngData->Init(theData);
141  hasFSData = true;
142  }
143  else if(dataType==12)
144  {
146  theFinalStatePhotons->InitMean(theData);
147  hasFSData = true;
148  }
149  else if(dataType==13)
150  {
153  hasFSData = true;
154  }
155  else if(dataType==14)
156  {
158  hasFSData = true;
159  }
160  else if(dataType==15)
161  {
163  hasFSData = true;
164  }
165  else
166  {
167  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4NeutronHPInelasticBaseFS");
168  }
169  }
170  //theData.close();
171 }
172 
174  G4ParticleDefinition ** theDefs,
175  G4int nDef)
176 {
177 
178 // prepare neutron
179  theResult.Clear();
180  G4double eKinetic = theTrack.GetKineticEnergy();
181  const G4HadProjectile *incidentParticle = &theTrack;
182  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
183  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
184  theNeutron.SetKineticEnergy( eKinetic );
185 
186 // prepare target
187  G4double targetMass;
188  G4double eps = 0.0001;
189  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
191 
192  if(theEnergyAngData!=0)
193  { targetMass = theEnergyAngData->GetTargetMass(); }
195  { targetMass = theAngularDistribution->GetTargetMass(); }
196 
197 //080731a
198 //110512 TKDB ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded.
199 //if ( targetMass == 0 ) G4cout << "080731a It looks like something wrong value in G4NDL, please update the latest version. If you use the latest, then please report this problem to Geant4 Hyper news." << G4endl;
200  if ( targetMass == 0 )
201  {
202  //G4cout << "TKDB targetMass = 0; ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded. This could be a similar situation." << G4endl;
203  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / G4Neutron::Neutron()->GetPDGMass();
204  }
205 
206  G4Nucleus aNucleus;
208  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
209  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
210 
211 // prepare energy in target rest frame
212  G4ReactionProduct boosted;
213  boosted.Lorentz(theNeutron, theTarget);
214  eKinetic = boosted.GetKineticEnergy();
215  G4double orgMomentum = boosted.GetMomentum().mag();
216 
217 // Take N-body phase-space distribution, if no other data present.
218  if(!HasFSData()) // adding the residual is trivial here @@@
219  {
220  G4NeutronHPNBodyPhaseSpace thePhaseSpaceDistribution;
221  G4double aPhaseMass=0;
222  G4int ii;
223  for(ii=0; ii<nDef; ii++)
224  {
225  aPhaseMass+=theDefs[ii]->GetPDGMass();
226  }
227  thePhaseSpaceDistribution.Init(aPhaseMass, nDef);
228  thePhaseSpaceDistribution.SetNeutron(&theNeutron);
229  thePhaseSpaceDistribution.SetTarget(&theTarget);
230  for(ii=0; ii<nDef; ii++)
231  {
232  G4double massCode = 1000.*std::abs(theDefs[ii]->GetPDGCharge());
233  massCode += theDefs[ii]->GetBaryonNumber();
234  G4double dummy = 0;
235  G4ReactionProduct * aSec = thePhaseSpaceDistribution.Sample(eKinetic, massCode, dummy);
236  aSec->Lorentz(*aSec, -1.*theTarget);
237  G4DynamicParticle * aPart = new G4DynamicParticle();
238  aPart->SetDefinition(aSec->GetDefinition());
239  aPart->SetMomentum(aSec->GetMomentum());
240  delete aSec;
241  theResult.AddSecondary(aPart);
242  }
244 
245  //TK120607
246  //Final momentum check should be done before return
248  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
249  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
250  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
251  adjust_final_state ( init_4p_lab );
252 
253  return;
254  }
255 
256 // set target and neutron in the relevant exit channel
257  if(theAngularDistribution!=0)
258  {
259  theAngularDistribution->SetTarget(theTarget);
260  theAngularDistribution->SetNeutron(theNeutron);
261  }
262  else if(theEnergyAngData!=0)
263  {
264  theEnergyAngData->SetTarget(theTarget);
265  theEnergyAngData->SetNeutron(theNeutron);
266  }
267 
268  G4ReactionProductVector * tmpHadrons = 0;
269  G4int ii, dummy;
270  unsigned int i;
271  if(theEnergyAngData != 0)
272  {
273  tmpHadrons = theEnergyAngData->Sample(eKinetic);
274  }
275  else if(theAngularDistribution!= 0)
276  {
277  G4bool * Done = new G4bool[nDef];
278  G4int i0;
279  for(i0=0; i0<nDef; i0++) Done[i0] = false;
280  if(tmpHadrons == 0)
281  {
282  tmpHadrons = new G4ReactionProductVector;
283  }
284  else
285  {
286  for(i=0; i<tmpHadrons->size(); i++)
287  {
288  for(ii=0; ii<nDef; ii++)
289  if(!Done[ii] && tmpHadrons->operator[](i)->GetDefinition() == theDefs[ii])
290  Done[ii] = true;
291  }
292  }
293  G4ReactionProduct * aHadron;
294  G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)));
295  G4ThreeVector bufferedDirection(0,0,0);
296  for(i0=0; i0<nDef; i0++)
297  {
298  if(!Done[i0])
299  {
300  aHadron = new G4ReactionProduct;
301  if(theEnergyDistribution!=0)
302  {
303  aHadron->SetDefinition(theDefs[i0]);
304  aHadron->SetKineticEnergy(theEnergyDistribution->Sample(eKinetic, dummy));
305  }
306  else if(nDef == 1)
307  {
308  aHadron->SetDefinition(theDefs[i0]);
309  aHadron->SetKineticEnergy(eKinetic);
310  }
311  else if(nDef == 2)
312  {
313  aHadron->SetDefinition(theDefs[i0]);
314  aHadron->SetKineticEnergy(50*MeV);
315  }
316  else
317  {
318  throw G4HadronicException(__FILE__, __LINE__, "No energy distribution to sample from in InelasticBaseFS::BaseApply");
319  }
321  if(theEnergyDistribution==0 && nDef == 2)
322  {
323  if(i0==0)
324  {
325  G4double mass1 = theDefs[0]->GetPDGMass();
326  G4double mass2 = theDefs[1]->GetPDGMass();
328  G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge());
329  G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber();
330  G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1);
331  G4double availableEnergy = eKinetic+massn+localMass-mass1-mass2-concreteMass;
332  // available kinetic energy in CMS (non relativistic)
333  G4double emin = availableEnergy+mass1+mass2 - std::sqrt((mass1+mass2)*(mass1+mass2)+orgMomentum*orgMomentum);
334  G4double p1=std::sqrt(2.*mass2*emin);
335  bufferedDirection = p1*aHadron->GetMomentum().unit();
336  if(getenv("HTOKEN")) // @@@@@ verify the nucleon counting...
337  {
338  G4cout << "HTOKEN "<<z1<<" "<<theBaseZ<<" "<<a1<<" "<<theBaseA<<" "<<availableEnergy<<" "
339  << emin<<G4endl;
340  }
341  }
342  else
343  {
344  bufferedDirection = -bufferedDirection;
345  }
346  // boost from cms to lab
347  if(getenv("HTOKEN"))
348  {
349  G4cout << " HTOKEN "<<bufferedDirection.mag2()<<G4endl;
350  }
351  aHadron->SetTotalEnergy( std::sqrt(aHadron->GetMass()*aHadron->GetMass()
352  +bufferedDirection.mag2()) );
353  aHadron->SetMomentum(bufferedDirection);
354  aHadron->Lorentz(*aHadron, -1.*(theTarget+theNeutron));
355  if(getenv("HTOKEN"))
356  {
357  G4cout << " HTOKEN "<<aHadron->GetTotalEnergy()<<" "<<aHadron->GetMomentum()<<G4endl;
358  }
359  }
360  tmpHadrons->push_back(aHadron);
361  }
362  }
363  delete [] Done;
364  }
365  else
366  {
367  throw G4HadronicException(__FILE__, __LINE__, "No data to create the neutrons in NInelasticFS");
368  }
369 
370  G4ReactionProductVector * thePhotons = 0;
371  if(theFinalStatePhotons!=0)
372  {
373  // the photon distributions are in the Nucleus rest frame.
374  G4ReactionProduct boosted_tmp;
375  boosted_tmp.Lorentz(theNeutron, theTarget);
376  G4double anEnergy = boosted_tmp.GetKineticEnergy();
377  thePhotons = theFinalStatePhotons->GetPhotons(anEnergy);
378  if(thePhotons!=0)
379  {
380  for(i=0; i<thePhotons->size(); i++)
381  {
382  // back to lab
383  thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), -1.*theTarget);
384  }
385  }
386  }
387  else if(theEnergyAngData!=0)
388  {
389 
390  G4double theGammaEnergy = theEnergyAngData->GetTotalMeanEnergy();
391  G4double anEnergy = boosted.GetKineticEnergy();
392  theGammaEnergy = anEnergy-theGammaEnergy;
393  theGammaEnergy += theNuclearMassDifference;
394  G4double eBindProducts = 0;
395  G4double eBindN = 0;
396  G4double eBindP = 0;
401  G4int ia=0;
402  for(i=0; i<tmpHadrons->size(); i++)
403  {
404  if(tmpHadrons->operator[](i)->GetDefinition() == G4Neutron::Neutron())
405  {
406  eBindProducts+=eBindN;
407  }
408  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Proton::Proton())
409  {
410  eBindProducts+=eBindP;
411  }
412  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Deuteron::Deuteron())
413  {
414  eBindProducts+=eBindD;
415  }
416  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Triton::Triton())
417  {
418  eBindProducts+=eBindT;
419  }
420  else if(tmpHadrons->operator[](i)->GetDefinition() == G4He3::He3())
421  {
422  eBindProducts+=eBindHe3;
423  }
424  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Alpha::Alpha())
425  {
426  eBindProducts+=eBindA;
427  ia++;
428  }
429  }
430 
431  theGammaEnergy += eBindProducts;
432 
433 //101111
434 //Special treatment for Be9 + n -> 2n + Be8 -> 2n + a + a
435 if ( (G4int)(theBaseZ+eps) == 4 && (G4int)(theBaseA+eps) == 9 )
436 {
437  // This only valid for G4NDL3.13,,,
438  if ( std::abs( theNuclearMassDifference -
440  G4NucleiProperties::GetBindingEnergy( 9 , 4 ) ) ) < 1*keV
441  && ia == 2 )
442  {
443  theGammaEnergy -= (2*eBindA);
444  }
445 }
446 
447  G4ReactionProductVector * theOtherPhotons = 0;
448  G4int iLevel;
449  while(theGammaEnergy>=theGammas.GetLevelEnergy(0))
450  {
451  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
452  {
453  if(theGammas.GetLevelEnergy(iLevel)<theGammaEnergy) break;
454  }
455  if(iLevel==0||iLevel==theGammas.GetNumberOfLevels()-1)
456  {
457  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
458  }
459  else
460  {
461  G4double random = G4UniformRand();
462  G4double eLow = theGammas.GetLevelEnergy(iLevel);
463  G4double eHigh = theGammas.GetLevelEnergy(iLevel+1);
464  if(random > (eHigh-eLow)/(theGammaEnergy-eLow)) iLevel++;
465  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
466  }
467  if(thePhotons==0) thePhotons = new G4ReactionProductVector;
468  if(theOtherPhotons != 0)
469  {
470  for(unsigned int iii=0; iii<theOtherPhotons->size(); iii++)
471  {
472  thePhotons->push_back(theOtherPhotons->operator[](iii));
473  }
474  delete theOtherPhotons;
475  }
476  theGammaEnergy -= theGammas.GetLevelEnergy(iLevel);
477  if(iLevel == -1) break;
478  }
479  }
480 
481 // fill the result
482  unsigned int nSecondaries = tmpHadrons->size();
483  unsigned int nPhotons = 0;
484  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
485  nSecondaries += nPhotons;
486  G4DynamicParticle * theSec;
487 
488  for(i=0; i<nSecondaries-nPhotons; i++)
489  {
490  theSec = new G4DynamicParticle;
491  theSec->SetDefinition(tmpHadrons->operator[](i)->GetDefinition());
492  theSec->SetMomentum(tmpHadrons->operator[](i)->GetMomentum());
493  theResult.AddSecondary(theSec);
494  delete tmpHadrons->operator[](i);
495  }
496  if(thePhotons != 0)
497  {
498  for(i=0; i<nPhotons; i++)
499  {
500  theSec = new G4DynamicParticle;
501  theSec->SetDefinition(thePhotons->operator[](i)->GetDefinition());
502  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
503  theResult.AddSecondary(theSec);
504  delete thePhotons->operator[](i);
505  }
506  }
507 
508 // some garbage collection
509  delete thePhotons;
510  delete tmpHadrons;
511 
512 //080721
514  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
515  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
516  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
517  adjust_final_state ( init_4p_lab );
518 
519 // clean up the primary neutron
521 }
G4ReactionProductVector * GetPhotons(G4double anEnergy)
void SetNeutron(const G4ReactionProduct &aNeutron)
void SetTarget(G4ReactionProduct &aTarget)
G4NeutronHPPhotonDist * theFinalStatePhotons
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetMomentum(const G4ThreeVector &momentum)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void Init(G4double aMass, G4int aCount)
void SetNeutron(G4ReactionProduct *aNeutron)
void InitGammas(G4double AR, G4double ZR)
void SetKineticEnergy(const G4double en)
static G4NeutronHPManager * GetInstance()
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
void Init(std::istream &aDataFile)
void SetTarget(const G4ReactionProduct &aTarget)
G4double Sample(G4double anEnergy, G4int &it)
void SetTarget(G4ReactionProduct *aTarget)
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
void GetDataStream(G4String, std::istringstream &iss)
G4NeutronHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
G4NeutronHPEnergyDistribution * theEnergyDistribution
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * GetDefinition() const
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void InitAngular(std::istream &aDataFile)
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
void SetNeutron(G4ReactionProduct &aNeutron)
G4double GetKineticEnergy() const
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &bit)
void SetTotalEnergy(const G4double en)
G4ErrorTarget * theTarget
Definition: errprop.cc:59
void InitEnergies(std::istream &aDataFile)
static G4Triton * Triton()
Definition: G4Triton.cc:95
void BaseApply(const G4HadProjectile &theTrack, G4ParticleDefinition **theDefs, G4int nDef)
static G4Proton * Proton()
Definition: G4Proton.cc:93
void Init(std::istream &aDataFile)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4NeutronHPEnAngCorrelation * theEnergyAngData
const G4LorentzVector & Get4Momentum() const
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:80
#define INT_MAX
Definition: templates.hh:111
G4double GetKineticEnergy() const
void SampleAndUpdate(G4ReactionProduct &aNeutron)
G4double GetTotalEnergy() const
G4bool InitMean(std::istream &aDataFile)
G4double total(Particle const *const p1, Particle const *const p2)
G4double GetPDGMass() const
static G4double GetBindingEnergy(const G4int A, const G4int Z)
G4double GetLevelEnergy(G4int aLevel)
Hep3Vector unit() const
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
double mag2() const
G4ThreeVector GetMomentum() const
G4NeutronHPAngular * theAngularDistribution
G4double GetTemperature() const
Definition: G4Material.hh:180
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
G4ReactionProductVector * Sample(G4double anEnergy)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
G4ReactionProductVector * GetDecayGammas(G4int aLevel)
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4NeutronHPDataUsed used)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void adjust_final_state(G4LorentzVector)
G4double GetPDGCharge() const
static G4He3 * He3()
Definition: G4He3.cc:94
double mag() const
G4double GetMass() const
void InitPartials(std::istream &aDataFile)
void AddSecondary(G4DynamicParticle *aP)
void Init(std::istream &aDataFile)