Geant4  9.6.p02
 All Classes 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 "G4SystemOfUnits.hh"
37 #include "G4Nucleus.hh"
38 #include "G4NucleiProperties.hh"
39 #include "G4He3.hh"
40 #include "G4Alpha.hh"
41 #include "G4Electron.hh"
42 #include "G4NeutronHPDataUsed.hh"
43 
44 #include "G4ParticleTable.hh"
45 
47 {
48  // char the[100] = {""};
49  // std::ostrstream ost(the, 100, std::ios::out);
50  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
51  // G4String * aName = new G4String(the);
52  // std::ifstream from(*aName, std::ios::in);
53 
54  std::ostringstream ost;
55  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
56  G4String aName = ost.str();
57  std::ifstream from(aName, std::ios::in);
58 
59  if(!from) return; // no data found for this isotope
60  // std::ifstream theGammaData(*aName, std::ios::in);
61  std::ifstream theGammaData(aName, std::ios::in);
62 
63  G4double eps = 0.001;
65  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) -
66  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
67  theGammas.Init(theGammaData);
68  // delete aName;
69 }
70 
72 {
73  gammaPath = "/Inelastic/Gammas/";
74  if(!getenv("G4NEUTRONHPDATA"))
75  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
76  G4String tBase = getenv("G4NEUTRONHPDATA");
77  gammaPath = tBase+gammaPath;
78  G4String tString = dirName;
79  G4bool dbool;
80  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M,tString, bit, dbool);
81  G4String filename = aFile.GetName();
82  SetAZMs( A, Z, M, aFile);
83  //theBaseA = aFile.GetA();
84  //theBaseZ = aFile.GetZ();
85  // theNDLDataA = (int)aFile.GetA();
86  // theNDLDataZ = aFile.GetZ();
87  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
88  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
89  {
90  if(getenv("NeutronHPNamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
91  hasAnyData = false;
92  hasFSData = false;
93  hasXsec = false;
94  return;
95  }
96  //theBaseA = A;
97  //theBaseZ = G4int(Z+.5);
98  std::ifstream theData(filename, std::ios::in);
99  if(!(theData))
100  {
101  hasAnyData = false;
102  hasFSData = false;
103  hasXsec = false;
104  theData.close();
105  return; // no data for exactly this isotope and FS
106  }
107  // here we go
108  G4int infoType, dataType, dummy=INT_MAX;
109  hasFSData = false;
110  while (theData >> infoType)
111  {
112  theData >> dataType;
113  if(dummy==INT_MAX) theData >> dummy >> dummy;
114  if(dataType==3)
115  {
116  G4int total;
117  theData >> total;
118  theXsection->Init(theData, total, eV);
119  }
120  else if(dataType==4)
121  {
123  theAngularDistribution->Init(theData);
124  hasFSData = true;
125  }
126  else if(dataType==5)
127  {
129  theEnergyDistribution->Init(theData);
130  hasFSData = true;
131  }
132  else if(dataType==6)
133  {
135  theEnergyAngData->Init(theData);
136  hasFSData = true;
137  }
138  else if(dataType==12)
139  {
141  theFinalStatePhotons->InitMean(theData);
142  hasFSData = true;
143  }
144  else if(dataType==13)
145  {
148  hasFSData = true;
149  }
150  else if(dataType==14)
151  {
153  hasFSData = true;
154  }
155  else if(dataType==15)
156  {
158  hasFSData = true;
159  }
160  else
161  {
162  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4NeutronHPInelasticBaseFS");
163  }
164  }
165  theData.close();
166 }
167 
169  G4ParticleDefinition ** theDefs,
170  G4int nDef)
171 {
172 
173 // prepare neutron
174  theResult.Clear();
175  G4double eKinetic = theTrack.GetKineticEnergy();
176  const G4HadProjectile *incidentParticle = &theTrack;
177  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
178  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
179  theNeutron.SetKineticEnergy( eKinetic );
180 
181 // prepare target
182  G4double targetMass;
183  G4double eps = 0.0001;
184  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
186 
187  if(theEnergyAngData!=0)
188  { targetMass = theEnergyAngData->GetTargetMass(); }
190  { targetMass = theAngularDistribution->GetTargetMass(); }
191 
192 //080731a
193 //110512 TKDB ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded.
194 //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;
195  if ( targetMass == 0 )
196  {
197  //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;
198  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / G4Neutron::Neutron()->GetPDGMass();
199  }
200 
201  G4Nucleus aNucleus;
203  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
204  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
205 
206 // prepare energy in target rest frame
207  G4ReactionProduct boosted;
208  boosted.Lorentz(theNeutron, theTarget);
209  eKinetic = boosted.GetKineticEnergy();
210  G4double orgMomentum = boosted.GetMomentum().mag();
211 
212 // Take N-body phase-space distribution, if no other data present.
213  if(!HasFSData()) // adding the residual is trivial here @@@
214  {
215  G4NeutronHPNBodyPhaseSpace thePhaseSpaceDistribution;
216  G4double aPhaseMass=0;
217  G4int ii;
218  for(ii=0; ii<nDef; ii++)
219  {
220  aPhaseMass+=theDefs[ii]->GetPDGMass();
221  }
222  thePhaseSpaceDistribution.Init(aPhaseMass, nDef);
223  thePhaseSpaceDistribution.SetNeutron(&theNeutron);
224  thePhaseSpaceDistribution.SetTarget(&theTarget);
225  for(ii=0; ii<nDef; ii++)
226  {
227  G4double massCode = 1000.*std::abs(theDefs[ii]->GetPDGCharge());
228  massCode += theDefs[ii]->GetBaryonNumber();
229  G4double dummy = 0;
230  G4ReactionProduct * aSec = thePhaseSpaceDistribution.Sample(eKinetic, massCode, dummy);
231  aSec->Lorentz(*aSec, -1.*theTarget);
232  G4DynamicParticle * aPart = new G4DynamicParticle();
233  aPart->SetDefinition(aSec->GetDefinition());
234  aPart->SetMomentum(aSec->GetMomentum());
235  delete aSec;
236  theResult.AddSecondary(aPart);
237  }
239 
240  //TK120607
241  //Final momentum check should be done before return
243  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
244  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
245  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
246  adjust_final_state ( init_4p_lab );
247 
248  return;
249  }
250 
251 // set target and neutron in the relevant exit channel
252  if(theAngularDistribution!=0)
253  {
254  theAngularDistribution->SetTarget(theTarget);
255  theAngularDistribution->SetNeutron(theNeutron);
256  }
257  else if(theEnergyAngData!=0)
258  {
259  theEnergyAngData->SetTarget(theTarget);
260  theEnergyAngData->SetNeutron(theNeutron);
261  }
262 
263  G4ReactionProductVector * tmpHadrons = 0;
264  G4int ii, dummy;
265  unsigned int i;
266  if(theEnergyAngData != 0)
267  {
268  tmpHadrons = theEnergyAngData->Sample(eKinetic);
269  }
270  else if(theAngularDistribution!= 0)
271  {
272  G4bool * Done = new G4bool[nDef];
273  G4int i0;
274  for(i0=0; i0<nDef; i0++) Done[i0] = false;
275  if(tmpHadrons == 0)
276  {
277  tmpHadrons = new G4ReactionProductVector;
278  }
279  else
280  {
281  for(i=0; i<tmpHadrons->size(); i++)
282  {
283  for(ii=0; ii<nDef; ii++)
284  if(!Done[ii] && tmpHadrons->operator[](i)->GetDefinition() == theDefs[ii])
285  Done[ii] = true;
286  }
287  }
288  G4ReactionProduct * aHadron;
289  G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)));
290  G4ThreeVector bufferedDirection(0,0,0);
291  for(i0=0; i0<nDef; i0++)
292  {
293  if(!Done[i0])
294  {
295  aHadron = new G4ReactionProduct;
296  if(theEnergyDistribution!=0)
297  {
298  aHadron->SetDefinition(theDefs[i0]);
299  aHadron->SetKineticEnergy(theEnergyDistribution->Sample(eKinetic, dummy));
300  }
301  else if(nDef == 1)
302  {
303  aHadron->SetDefinition(theDefs[i0]);
304  aHadron->SetKineticEnergy(eKinetic);
305  }
306  else if(nDef == 2)
307  {
308  aHadron->SetDefinition(theDefs[i0]);
309  aHadron->SetKineticEnergy(50*MeV);
310  }
311  else
312  {
313  throw G4HadronicException(__FILE__, __LINE__, "No energy distribution to sample from in InelasticBaseFS::BaseApply");
314  }
316  if(theEnergyDistribution==0 && nDef == 2)
317  {
318  if(i0==0)
319  {
320  G4double mass1 = theDefs[0]->GetPDGMass();
321  G4double mass2 = theDefs[1]->GetPDGMass();
323  G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge());
324  G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber();
325  G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1);
326  G4double availableEnergy = eKinetic+massn+localMass-mass1-mass2-concreteMass;
327  // available kinetic energy in CMS (non relativistic)
328  G4double emin = availableEnergy+mass1+mass2 - std::sqrt((mass1+mass2)*(mass1+mass2)+orgMomentum*orgMomentum);
329  G4double p1=std::sqrt(2.*mass2*emin);
330  bufferedDirection = p1*aHadron->GetMomentum().unit();
331  if(getenv("HTOKEN")) // @@@@@ verify the nucleon counting...
332  {
333  G4cout << "HTOKEN "<<z1<<" "<<theBaseZ<<" "<<a1<<" "<<theBaseA<<" "<<availableEnergy<<" "
334  << emin<<G4endl;
335  }
336  }
337  else
338  {
339  bufferedDirection = -bufferedDirection;
340  }
341  // boost from cms to lab
342  if(getenv("HTOKEN"))
343  {
344  G4cout << " HTOKEN "<<bufferedDirection.mag2()<<G4endl;
345  }
346  aHadron->SetTotalEnergy( std::sqrt(aHadron->GetMass()*aHadron->GetMass()
347  +bufferedDirection.mag2()) );
348  aHadron->SetMomentum(bufferedDirection);
349  aHadron->Lorentz(*aHadron, -1.*(theTarget+theNeutron));
350  if(getenv("HTOKEN"))
351  {
352  G4cout << " HTOKEN "<<aHadron->GetTotalEnergy()<<" "<<aHadron->GetMomentum()<<G4endl;
353  }
354  }
355  tmpHadrons->push_back(aHadron);
356  }
357  }
358  delete [] Done;
359  }
360  else
361  {
362  throw G4HadronicException(__FILE__, __LINE__, "No data to create the neutrons in NInelasticFS");
363  }
364 
365  G4ReactionProductVector * thePhotons = 0;
366  if(theFinalStatePhotons!=0)
367  {
368  // the photon distributions are in the Nucleus rest frame.
369  G4ReactionProduct boosted_tmp;
370  boosted_tmp.Lorentz(theNeutron, theTarget);
371  G4double anEnergy = boosted_tmp.GetKineticEnergy();
372  thePhotons = theFinalStatePhotons->GetPhotons(anEnergy);
373  if(thePhotons!=0)
374  {
375  for(i=0; i<thePhotons->size(); i++)
376  {
377  // back to lab
378  thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), -1.*theTarget);
379  }
380  }
381  }
382  else if(theEnergyAngData!=0)
383  {
384 
385  G4double theGammaEnergy = theEnergyAngData->GetTotalMeanEnergy();
386  G4double anEnergy = boosted.GetKineticEnergy();
387  theGammaEnergy = anEnergy-theGammaEnergy;
388  theGammaEnergy += theNuclearMassDifference;
389  G4double eBindProducts = 0;
390  G4double eBindN = 0;
391  G4double eBindP = 0;
396  G4int ia=0;
397  for(i=0; i<tmpHadrons->size(); i++)
398  {
399  if(tmpHadrons->operator[](i)->GetDefinition() == G4Neutron::Neutron())
400  {
401  eBindProducts+=eBindN;
402  }
403  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Proton::Proton())
404  {
405  eBindProducts+=eBindP;
406  }
407  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Deuteron::Deuteron())
408  {
409  eBindProducts+=eBindD;
410  }
411  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Triton::Triton())
412  {
413  eBindProducts+=eBindT;
414  }
415  else if(tmpHadrons->operator[](i)->GetDefinition() == G4He3::He3())
416  {
417  eBindProducts+=eBindHe3;
418  }
419  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Alpha::Alpha())
420  {
421  eBindProducts+=eBindA;
422  ia++;
423  }
424  }
425 
426  theGammaEnergy += eBindProducts;
427 
428 //101111
429 //Special treatment for Be9 + n -> 2n + Be8 -> 2n + a + a
430 if ( (G4int)(theBaseZ+eps) == 4 && (G4int)(theBaseA+eps) == 9 )
431 {
432  // This only valid for G4NDL3.13,,,
433  if ( std::abs( theNuclearMassDifference -
435  G4NucleiProperties::GetBindingEnergy( 9 , 4 ) ) ) < 1*keV
436  && ia == 2 )
437  {
438  theGammaEnergy -= (2*eBindA);
439  }
440 }
441 
442  G4ReactionProductVector * theOtherPhotons = 0;
443  G4int iLevel;
444  while(theGammaEnergy>=theGammas.GetLevelEnergy(0))
445  {
446  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
447  {
448  if(theGammas.GetLevelEnergy(iLevel)<theGammaEnergy) break;
449  }
450  if(iLevel==0||iLevel==theGammas.GetNumberOfLevels()-1)
451  {
452  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
453  }
454  else
455  {
456  G4double random = G4UniformRand();
457  G4double eLow = theGammas.GetLevelEnergy(iLevel);
458  G4double eHigh = theGammas.GetLevelEnergy(iLevel+1);
459  if(random > (eHigh-eLow)/(theGammaEnergy-eLow)) iLevel++;
460  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
461  }
462  if(thePhotons==0) thePhotons = new G4ReactionProductVector;
463  if(theOtherPhotons != 0)
464  {
465  for(unsigned int iii=0; iii<theOtherPhotons->size(); iii++)
466  {
467  thePhotons->push_back(theOtherPhotons->operator[](iii));
468  }
469  delete theOtherPhotons;
470  }
471  theGammaEnergy -= theGammas.GetLevelEnergy(iLevel);
472  if(iLevel == -1) break;
473  }
474  }
475 
476 // fill the result
477  unsigned int nSecondaries = tmpHadrons->size();
478  unsigned int nPhotons = 0;
479  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
480  nSecondaries += nPhotons;
481  G4DynamicParticle * theSec;
482 
483  for(i=0; i<nSecondaries-nPhotons; i++)
484  {
485  theSec = new G4DynamicParticle;
486  theSec->SetDefinition(tmpHadrons->operator[](i)->GetDefinition());
487  theSec->SetMomentum(tmpHadrons->operator[](i)->GetMomentum());
488  theResult.AddSecondary(theSec);
489  delete tmpHadrons->operator[](i);
490  }
491  if(thePhotons != 0)
492  {
493  for(i=0; i<nPhotons; i++)
494  {
495  theSec = new G4DynamicParticle;
496  theSec->SetDefinition(thePhotons->operator[](i)->GetDefinition());
497  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
498  theResult.AddSecondary(theSec);
499  delete thePhotons->operator[](i);
500  }
501  }
502 
503 // some garbage collection
504  delete thePhotons;
505  delete tmpHadrons;
506 
507 //080721
509  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
510  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
511  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
512  adjust_final_state ( init_4p_lab );
513 
514 // clean up the primary neutron
516 }