ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParticleHPInelasticBaseFS.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4ParticleHPInelasticBaseFS.cc
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 // particle_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 //
35 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
36 //
37 // June-2019 - E. Mendoza --> Added protection against residual with Z<0 or A<Z + adjust_final_state is not applied when data is in MF=6 format (no correlated particle emission) + bug correction (add Q value info to G4ParticleHPNBodyPhaseSpace).
38 
39 
41 #include "G4ParticleHPManager.hh"
42 #include "G4Nucleus.hh"
43 #include "G4NucleiProperties.hh"
44 #include "G4He3.hh"
45 #include "G4Alpha.hh"
46 #include "G4Electron.hh"
47 #include "G4ParticleHPDataUsed.hh"
48 
49 #include "G4IonTable.hh"
50 
52 {
53  // char the[100] = {""};
54  // std::ostrstream ost(the, 100, std::ios::out);
55  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
56  // G4String * aName = new G4String(the);
57  // std::ifstream from(*aName, std::ios::in);
58 
59  std::ostringstream ost;
60  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
61  G4String aName = ost.str();
62  std::ifstream from(aName, std::ios::in);
63 
64  if(!from) return; // no data found for this isotope
65  // std::ifstream theGammaData(*aName, std::ios::in);
66  std::ifstream theGammaData(aName, std::ios::in);
67 
68  G4double eps = 0.001;
70  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) -
71  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
72  theGammas.Init(theGammaData);
73  // delete aName;
74 
75 }
76 
78 {
79  gammaPath = "/Inelastic/Gammas/";
80  if(!std::getenv("G4NEUTRONHPDATA"))
81  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files where Inelastic/Gammas data is found.");
82  G4String tBase = std::getenv("G4NEUTRONHPDATA");
83  gammaPath = tBase+gammaPath;
84  G4String tString = dirName;
85  G4bool dbool;
86  G4ParticleHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M,tString, bit, dbool);
87  G4String filename = aFile.GetName();
88 #ifdef G4PHPDEBUG
89  if( std::getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPInelasticBaseFS::Init FILE " << filename << G4endl;
90 #endif
91  SetAZMs( A, Z, M, aFile);
92  //theBaseA = aFile.GetA();
93  //theBaseZ = aFile.GetZ();
94  // theNDLDataA = (int)aFile.GetA();
95  // theNDLDataZ = aFile.GetZ();
96  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
97  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
98  {
99 #ifdef G4PHPDEBUG
100  if(std::getenv("G4ParticleHPDebug_NamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
101 #endif
102  hasAnyData = false;
103  hasFSData = false;
104  hasXsec = false;
105  return;
106  }
107  // theBaseA = A;
108  // theBaseZ = G4int(Z+.5);
109  //std::ifstream theData(filename, std::ios::in);
110  //130205 For compressed data files
111  std::istringstream theData(std::ios::in);
112  G4ParticleHPManager::GetInstance()->GetDataStream(filename,theData);
113  //130205 END
114  if(!theData) //"!" is a operator of ios
115  {
116  hasAnyData = false;
117  hasFSData = false;
118  hasXsec = false;
119  // theData.close();
120  return; // no data for exactly this isotope and FS
121  }
122  // here we go
123  G4int infoType, dataType, dummy=INT_MAX;
124  hasFSData = false;
125  while (theData >> infoType) // Loop checking, 11.05.2015, T. Koi
126  {
127  theData >> dataType;
128 
129  if(dummy==INT_MAX) theData >> Qvalue >> dummy;
130  Qvalue*=CLHEP::eV; //In G4NDL4.5 this value is the MT number (<1000), in others is que Q-value in eV
131 
132  if(dataType==3)
133  {
134  G4int total;
135  theData >> total;
136  theXsection->Init(theData, total, CLHEP::eV);
137  }
138  else if(dataType==4)
139  {
141  theAngularDistribution->Init(theData);
142  hasFSData = true;
143  }
144  else if(dataType==5)
145  {
147  theEnergyDistribution->Init(theData);
148  hasFSData = true;
149  }
150  else if(dataType==6)
151  {
153  //- G4cout << this << " BaseFS theEnergyAngData " << theEnergyAngData << G4endl; //GDEB
154  theEnergyAngData->Init(theData);
155  hasFSData = true;
156  }
157  else if(dataType==12)
158  {
160  theFinalStatePhotons->InitMean(theData);
161  hasFSData = true;
162  }
163  else if(dataType==13)
164  {
167  hasFSData = true;
168  }
169  else if(dataType==14)
170  {
172  hasFSData = true;
173  }
174  else if(dataType==15)
175  {
177  hasFSData = true;
178  }
179  else
180  {
181  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4ParticleHPInelasticBaseFS");
182  }
183  }
184  // theData.close();
185 }
186 
188  G4ParticleDefinition ** theDefs,
189  G4int nDef)
190 {
191 
192 // prepare neutron
193  if ( theResult.Get() == NULL ) theResult.Put( new G4HadFinalState );
194  theResult.Get()->Clear();
195  G4double eKinetic = theTrack.GetKineticEnergy();
196  const G4HadProjectile *hadProjectile = &theTrack;
197  G4ReactionProduct incidReactionProduct( const_cast<G4ParticleDefinition *>(hadProjectile->GetDefinition()) );
198  incidReactionProduct.SetMomentum( hadProjectile->Get4Momentum().vect() );
199  incidReactionProduct.SetKineticEnergy( eKinetic );
200 
201 // prepare target
202  G4double targetMass;
203  G4double eps = 0.0001;
204  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
205  //theProjectile->GetPDGMass();
207 
208  //give priority to ENDF vales for target mass
209  if(theEnergyAngData!=0)
210  { targetMass = theEnergyAngData->GetTargetMass(); }
212  { targetMass = theAngularDistribution->GetTargetMass(); }
213 
214 //080731a
215 //110512 TKDB ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded.
216 //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;
217  if ( targetMass == 0 )
218  {
219  //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;
220  //targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / theProjectile->GetPDGMass();
221  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / G4Neutron::Neutron()->GetPDGMass();
222  }
223 
224  G4Nucleus aNucleus;
226  //G4ThreeVector neuVelo = (1./hadProjectile->GetDefinition()->GetPDGMass())*incidReactionProduct.GetMomentum();
227  //theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
228  //G4Nucleus::GetBiasedThermalNucleus requests normalization of mass and velocity in neutron mass
229  G4ThreeVector neuVelo = (1./G4Neutron::Neutron()->GetPDGMass())*incidReactionProduct.GetMomentum();
230  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
231 
232  theTarget.SetDefinition( G4IonTable::GetIonTable()->GetIon( G4int(theBaseZ), G4int(theBaseA) , 0.0 ) );
233 
234 // prepare energy in target rest frame
235  G4ReactionProduct boosted;
236  boosted.Lorentz(incidReactionProduct, theTarget);
237  eKinetic = boosted.GetKineticEnergy();
238  G4double orgMomentum = boosted.GetMomentum().mag();
239 
240 // Take N-body phase-space distribution, if no other data present.
241  if(!HasFSData()) // adding the residual is trivial here @@@
242  {
243  G4ParticleHPNBodyPhaseSpace thePhaseSpaceDistribution;
244  G4double aPhaseMass=0;
245  G4int ii;
246  for(ii=0; ii<nDef; ii++)
247  {
248  aPhaseMass+=theDefs[ii]->GetPDGMass();
249  }
250 
251  //----------------------------------------------------------------------------
252  if(Qvalue<1.*CLHEP::keV && Qvalue>-1.*CLHEP::keV){ //Not in the G4NDL lib or not calculated yet:
253  //Calculate residual:
254  G4int ResidualA=theBaseA;
255  G4int ResidualZ=theBaseZ;
256  for (ii = 0; ii < nDef; ii++) {
257  ResidualZ -= theDefs[ii]->GetAtomicNumber();
258  ResidualA -= theDefs[ii]->GetBaryonNumber();
259  }
260 
261  if (ResidualA > 0 && ResidualZ > 0) {
262  G4ParticleDefinition* resid = G4IonTable::GetIonTable()->GetIon(ResidualZ,ResidualA);
263  Qvalue = incidReactionProduct.GetMass()+theTarget.GetMass()-aPhaseMass-resid->GetPDGMass();
264  }
265 
266  if (Qvalue > 400*CLHEP::MeV || Qvalue < -400*CLHEP::MeV) {
267  //Then Q value is probably too large ...
268  Qvalue = 1.1*CLHEP::keV;
269  }
270  }
271  //----------------------------------------------------------------------------
272 
273  thePhaseSpaceDistribution.Init(aPhaseMass, nDef);
274  thePhaseSpaceDistribution.SetProjectileRP(&incidReactionProduct);
275  thePhaseSpaceDistribution.SetTarget(&theTarget);
276  thePhaseSpaceDistribution.SetQValue(Qvalue);
277 
278  for(ii=0; ii<nDef; ii++)
279  {
280  G4double massCode = 1000.*std::abs(theDefs[ii]->GetPDGCharge());
281  massCode += theDefs[ii]->GetBaryonNumber();
282  G4double dummy = 0;
283  G4ReactionProduct * aSec = thePhaseSpaceDistribution.Sample(eKinetic, massCode, dummy);
284  aSec->Lorentz(*aSec, -1.*theTarget);
285  G4DynamicParticle * aPart = new G4DynamicParticle();
286  aPart->SetDefinition(aSec->GetDefinition());
287  aPart->SetMomentum(aSec->GetMomentum());
288  delete aSec;
289  theResult.Get()->AddSecondary(aPart);
290 #ifdef G4PHPDEBUG
291  if( std::getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply NoFSData add secondary " << aPart->GetParticleDefinition()->GetParticleName() << " E= " << aPart->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
292 #endif
293  }
295  //TK120607
296  //Final momentum check should be done before return
298  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
299  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
300  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
301  adjust_final_state ( init_4p_lab );
302 
303  return;
304  }
305 
306 // set target and neutron in the relevant exit channel
307  if(theAngularDistribution!=0)
308  {
309  theAngularDistribution->SetTarget(theTarget);
310  theAngularDistribution->SetProjectileRP(incidReactionProduct);
311  }
312  else if(theEnergyAngData!=0)
313  {
314  theEnergyAngData->SetTarget(theTarget);
315  theEnergyAngData->SetProjectileRP(incidReactionProduct);
316  }
317 
318  G4ReactionProductVector * tmpHadrons = 0;
319 #ifdef G4PHPDEBUG
320  //To avoid compilation error around line 532.
321  G4int ii(0);
322 #endif
323  G4int dummy;
324  unsigned int i;
325 
326  if(theEnergyAngData != 0)
327  {
328  tmpHadrons = theEnergyAngData->Sample(eKinetic);
329 
330  if ( !std::getenv( "G4PHP_DO_NOT_ADJUST_FINAL_STATE" ) ) {
331  //141017 Fix BEGIN
332  //Adjust A and Z in the case of miss much between selected data and target nucleus
333  if ( tmpHadrons != NULL ) {
334  G4int sumA = 0;
335  G4int sumZ = 0;
336  G4int maxA = 0;
337  G4int jAtMaxA = 0;
338  for ( G4int j = 0 ; j != (G4int)tmpHadrons->size() ; j++ ) {
339  //G4cout << __FILE__ << " " << __LINE__ << "th line: tmpHadrons->at(j)->GetDefinition()->GetParticleName() = " << tmpHadrons->at(j)->GetDefinition()->GetParticleName() << G4endl;
340  if ( tmpHadrons->at(j)->GetDefinition()->GetBaryonNumber() > maxA ) {
341  maxA = tmpHadrons->at(j)->GetDefinition()->GetBaryonNumber();
342  jAtMaxA = j;
343  }
344  sumA += tmpHadrons->at(j)->GetDefinition()->GetBaryonNumber();
345  sumZ += G4int( tmpHadrons->at(j)->GetDefinition()->GetPDGCharge() + eps );
346  }
347  G4int dA = (G4int)theBaseA + hadProjectile->GetDefinition()->GetBaryonNumber() - sumA;
348  G4int dZ = (G4int)theBaseZ + G4int( hadProjectile->GetDefinition()->GetPDGCharge() + eps ) - sumZ;
349  if ( dA < 0 || dZ < 0 ) {
350  G4int newA = tmpHadrons->at(jAtMaxA)->GetDefinition()->GetBaryonNumber() + dA ;
351  G4int newZ = G4int( tmpHadrons->at(jAtMaxA)->GetDefinition()->GetPDGCharge() + eps ) + dZ;
352  if(newA>newZ && newZ>0){
353  G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon ( newZ , newA );
354  tmpHadrons->at( jAtMaxA )->SetDefinition( pd );
355  }
356  }
357  }
358  //141017 Fix END
359  }
360  }
361  else if(theAngularDistribution!= 0)
362  {
363  G4bool * Done = new G4bool[nDef];
364  G4int i0;
365  for(i0=0; i0<nDef; i0++) Done[i0] = false;
366 
367 //Following lines are commented out to fix coverity defeat 58622
368 // if(tmpHadrons == 0)
369 // {
370  tmpHadrons = new G4ReactionProductVector;
371 // }
372 // else
373 // {
374 // for(i=0; i<tmpHadrons->size(); i++)
375 // {
376 // for(ii=0; ii<nDef; ii++)
377 // if(!Done[ii] && tmpHadrons->operator[](i)->GetDefinition() == theDefs[ii])
378 // Done[ii] = true;
379 // }
380 #ifdef G4PHPDEBUG
381 // if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply secondary previously added " << tmpHadrons->operator[](i)->GetDefinition()->GetParticleName() << " E= " << tmpHadrons->operator[](i)->GetKineticEnergy() << G4endl;
382 #endif
383 // }
384  G4ReactionProduct * aHadron;
385  G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)));
386  G4ThreeVector bufferedDirection(0,0,0);
387  for(i0=0; i0<nDef; i0++)
388  {
389  if(!Done[i0])
390  {
391  aHadron = new G4ReactionProduct;
392  if(theEnergyDistribution!=0)
393  {
394  aHadron->SetDefinition(theDefs[i0]);
395  aHadron->SetKineticEnergy(theEnergyDistribution->Sample(eKinetic, dummy));
396  }
397  else if(nDef == 1)
398  {
399  aHadron->SetDefinition(theDefs[i0]);
400  aHadron->SetKineticEnergy(eKinetic);
401  }
402  else if(nDef == 2)
403  {
404  aHadron->SetDefinition(theDefs[i0]);
405  aHadron->SetKineticEnergy(50*CLHEP::MeV);
406  }
407  else
408  {
409  throw G4HadronicException(__FILE__, __LINE__, "No energy distribution to sample from in InelasticBaseFS::BaseApply");
410  }
412  if(theEnergyDistribution==0 && nDef == 2)
413  {
414  if(i0==0)
415  {
416  G4double mass1 = theDefs[0]->GetPDGMass();
417  G4double mass2 = theDefs[1]->GetPDGMass();
418  G4double massn = theProjectile->GetPDGMass();
419  G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge());
420  G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber();
421  G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1);
422  G4double availableEnergy = eKinetic+massn+localMass-mass1-mass2-concreteMass;
423  // available kinetic energy in CMS (non relativistic)
424  G4double emin = availableEnergy+mass1+mass2 - std::sqrt((mass1+mass2)*(mass1+mass2)+orgMomentum*orgMomentum);
425  G4double p1=std::sqrt(2.*mass2*emin);
426  bufferedDirection = p1*aHadron->GetMomentum().unit();
427 #ifdef G4PHPDEBUG
428  if(std::getenv("G4ParticleHPDebug")) // @@@@@ verify the nucleon counting...
429  {
430  G4cout << "G4ParticleHPInelasticBaseFS "<<z1<<" "<<theBaseZ<<" "<<a1<<" "<<theBaseA<<" "<<availableEnergy<<" "
431  << emin<<G4endl;
432  }
433 #endif
434  }
435  else
436  {
437  bufferedDirection = -bufferedDirection;
438  }
439  // boost from cms to lab
440 #ifdef G4PHPDEBUG
441  if(std::getenv("G4ParticleHPDebug"))
442  {
443  G4cout << " G4ParticleHPInelasticBaseFS "<<bufferedDirection.mag2()<<G4endl;
444  }
445 #endif
446  aHadron->SetTotalEnergy( std::sqrt(aHadron->GetMass()*aHadron->GetMass()
447  +bufferedDirection.mag2()) );
448  aHadron->SetMomentum(bufferedDirection);
449  aHadron->Lorentz(*aHadron, -1.*(theTarget+incidReactionProduct));
450 #ifdef G4PHPDEBUG
451  if(std::getenv("G4ParticleHPDebug"))
452  {
453  G4cout << " G4ParticleHPInelasticBaseFS "<<aHadron->GetTotalEnergy()<<" "<<aHadron->GetMomentum()<<G4endl;
454  }
455 #endif
456  }
457  tmpHadrons->push_back(aHadron);
458 #ifdef G4PHPDEBUG
459  if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply FSData add secondary " << aHadron->GetDefinition()->GetParticleName() << " E= " << aHadron->GetKineticEnergy() << G4endl;
460 #endif
461  }
462  }
463  delete [] Done;
464  }
465  else
466  {
467  throw G4HadronicException(__FILE__, __LINE__, "No data to create the neutrons in NInelasticFS");
468  }
469 
470  G4ReactionProductVector * thePhotons = 0;
471  if(theFinalStatePhotons!=0)
472  {
473  // the photon distributions are in the Nucleus rest frame.
474  G4ReactionProduct boosted_tmp;
475  boosted_tmp.Lorentz(incidReactionProduct, theTarget);
476  G4double anEnergy = boosted_tmp.GetKineticEnergy();
477  thePhotons = theFinalStatePhotons->GetPhotons(anEnergy);
478  if(thePhotons!=0)
479  {
480  for(i=0; i<thePhotons->size(); i++)
481  {
482  // back to lab
483  thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), -1.*theTarget);
484  }
485  }
486  }
487  else if(theEnergyAngData!=0)
488  {
489 
490  // PA130927: do not create photons to adjust binding energy
491  G4bool bAdjustPhotons = true;
492 #ifdef PHP_AS_HP
493  bAdjustPhotons = true;
494 #else
495  if ( std::getenv( "G4PHP_DO_NOT_ADJUST_FINAL_STATE" ) ) bAdjustPhotons = false;
496 #endif
497 
498  if( bAdjustPhotons ) {
499  G4double theGammaEnergy = theEnergyAngData->GetTotalMeanEnergy();
500  G4double anEnergy = boosted.GetKineticEnergy();
501  theGammaEnergy = anEnergy-theGammaEnergy;
502  theGammaEnergy += theNuclearMassDifference;
503  G4double eBindProducts = 0;
504  G4double eBindN = 0;
505  G4double eBindP = 0;
510  G4int ia=0;
511  for(i=0; i<tmpHadrons->size(); i++)
512  {
513  if(tmpHadrons->operator[](i)->GetDefinition() == G4Neutron::Neutron())
514  {
515  eBindProducts+=eBindN;
516  }
517  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Proton::Proton())
518  {
519  eBindProducts+=eBindP;
520  }
521  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Deuteron::Deuteron())
522  {
523  eBindProducts+=eBindD;
524  }
525  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Triton::Triton())
526  {
527  eBindProducts+=eBindT;
528  }
529  else if(tmpHadrons->operator[](i)->GetDefinition() == G4He3::He3())
530  {
531  eBindProducts+=eBindHe3;
532  }
533  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Alpha::Alpha())
534  {
535  eBindProducts+=eBindA;
536  ia++;
537  }
538  }
539 
540  theGammaEnergy += eBindProducts;
541 
542 #ifdef G4PHPDEBUG
543  if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply gamma Energy " << theGammaEnergy << " eBindProducts " << eBindProducts << G4endl;
544 #endif
545 
546  //101111
547  //Special treatment for Be9 + n -> 2n + Be8 -> 2n + a + a
548  if ( (G4int)(theBaseZ+eps) == 4 && (G4int)(theBaseA+eps) == 9 )
549  {
550  // This only valid for G4NDL3.13,,,
554  && ia == 2 )
555  {
556  theGammaEnergy -= (2*eBindA);
557  }
558  }
559 
560  G4ReactionProductVector * theOtherPhotons = 0;
561  G4int iLevel;
562  while(theGammaEnergy>=theGammas.GetLevelEnergy(0)) // Loop checking, 11.05.2015, T. Koi
563  {
564  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
565  {
566  if(theGammas.GetLevelEnergy(iLevel)<theGammaEnergy) break;
567  }
568  if(iLevel==0||iLevel==theGammas.GetNumberOfLevels()-1)
569  {
570  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
571 #ifdef G4PHPDEBUG
572  if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply adding gamma from level " << iLevel << theOtherPhotons->operator[](ii)->GetKineticEnergy() << G4endl;
573 #endif
574  }
575  else
576  {
577  G4double random = G4UniformRand();
578  G4double eLow = theGammas.GetLevelEnergy(iLevel);
579  G4double eHigh = theGammas.GetLevelEnergy(iLevel+1);
580  if(random > (eHigh-eLow)/(theGammaEnergy-eLow)) iLevel++;
581  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
582  }
583  if(thePhotons==0) thePhotons = new G4ReactionProductVector;
584  if(theOtherPhotons != 0)
585  {
586  for(unsigned int iii=0; iii<theOtherPhotons->size(); iii++)
587  {
588  thePhotons->push_back(theOtherPhotons->operator[](iii));
589 #ifdef G4PHPDEBUG
590  if( std::getenv("G4ParticleHPDebug"))
591  G4cout << iii << " G4ParticleHPInelasticBaseFS::BaseApply adding gamma " << theOtherPhotons->operator[](iii)->GetKineticEnergy() << G4endl;
592 #endif
593  }
594  delete theOtherPhotons;
595  }
596  theGammaEnergy -= theGammas.GetLevelEnergy(iLevel);
597  if(iLevel == -1) break;
598  }
599  }
600  }
601 
602  // fill the result
603  unsigned int nSecondaries = tmpHadrons->size();
604  unsigned int nPhotons = 0;
605  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
606  nSecondaries += nPhotons;
607  G4DynamicParticle * theSec;
608 #ifdef G4PHPDEBUG
609  if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply N hadrons " << nSecondaries-nPhotons << G4endl;
610 #endif
611 
612  for(i=0; i<nSecondaries-nPhotons; i++)
613  {
614  theSec = new G4DynamicParticle;
615  theSec->SetDefinition(tmpHadrons->operator[](i)->GetDefinition());
616  theSec->SetMomentum(tmpHadrons->operator[](i)->GetMomentum());
617  theResult.Get()->AddSecondary(theSec);
618 #ifdef G4PHPDEBUG
619  if( std::getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply add secondary2 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
620 #endif
621  if( std::getenv("G4PHPTEST") ) G4cout << " InelasticBaseFS COS THETA " << std::cos(theSec->GetMomentum().theta()) << " " << (theSec->GetMomentum().theta()) << " " << theSec->GetMomentum() << " E "<< theSec->GetKineticEnergy() << " " << theSec->GetDefinition()->GetParticleName() << G4endl; //GDEB
622  delete tmpHadrons->operator[](i);
623  }
624 #ifdef G4PHPDEBUG
625  if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply N photons " << nPhotons << G4endl;
626 #endif
627  if(thePhotons != 0)
628  {
629  for(i=0; i<nPhotons; i++)
630  {
631  theSec = new G4DynamicParticle;
632  theSec->SetDefinition(thePhotons->operator[](i)->GetDefinition());
633  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
634  theResult.Get()->AddSecondary(theSec);
635 #ifdef G4PHPDEBUG
636  if( std::getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply add secondary3 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
637 #endif
638  delete thePhotons->operator[](i);
639  }
640  }
641 
642 // some garbage collection
643  delete thePhotons;
644  delete tmpHadrons;
645 
646 //080721
648  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
649  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
650  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
651 
652  //if data in MF=6 format (no correlated particle emission), then adjust_final_state can give severe errors:
653  if(theEnergyAngData==0){adjust_final_state ( init_4p_lab );}
654 
655 // clean up the primary neutron
657 }