Geant4  10.01.p02
G4ParticleHPInelasticBaseFS.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 // 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 //
38 #include "G4ParticleHPManager.hh"
39 #include "G4Nucleus.hh"
40 #include "G4NucleiProperties.hh"
41 #include "G4He3.hh"
42 #include "G4Alpha.hh"
43 #include "G4Electron.hh"
44 #include "G4ParticleHPDataUsed.hh"
45 
46 #include "G4IonTable.hh"
47 
49 {
50  // char the[100] = {""};
51  // std::ostrstream ost(the, 100, std::ios::out);
52  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
53  // G4String * aName = new G4String(the);
54  // std::ifstream from(*aName, std::ios::in);
55 
56  std::ostringstream ost;
57  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
58  G4String aName = ost.str();
59  std::ifstream from(aName, std::ios::in);
60 
61  if(!from) return; // no data found for this isotope
62  // std::ifstream theGammaData(*aName, std::ios::in);
63  std::ifstream theGammaData(aName, std::ios::in);
64 
65  G4double eps = 0.001;
67  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) -
68  G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
69  theGammas.Init(theGammaData);
70  // delete aName;
71 
72 }
73 
75 {
76  gammaPath = "/Inelastic/Gammas/";
77  if(!getenv("G4NEUTRONHPDATA"))
78  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files where Inelastic/Gammas data is found.");
79  G4String tBase = getenv("G4NEUTRONHPDATA");
80  gammaPath = tBase+gammaPath;
81  G4String tString = dirName;
82  G4bool dbool;
83  G4ParticleHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M,tString, bit, dbool);
84  G4String filename = aFile.GetName();
85 #ifdef G4PHPDEBUG
86  if( getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPInelasticBaseFS::Init FILE " << filename << G4endl;
87 #endif
88  SetAZMs( A, Z, M, aFile);
89  //theBaseA = aFile.GetA();
90  //theBaseZ = aFile.GetZ();
91  // theNDLDataA = (int)aFile.GetA();
92  // theNDLDataZ = aFile.GetZ();
93  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
94  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
95  {
96 #ifdef G4PHPDEBUG
97  if(getenv("G4ParticleHPDebug_NamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
98 #endif
99  hasAnyData = false;
100  hasFSData = false;
101  hasXsec = false;
102  return;
103  }
104  // theBaseA = A;
105  // theBaseZ = G4int(Z+.5);
106  //std::ifstream theData(filename, std::ios::in);
107  //130205 For compressed data files
108  std::istringstream theData(std::ios::in);
109  G4ParticleHPManager::GetInstance()->GetDataStream(filename,theData);
110  //130205 END
111  if(!theData) //"!" is a operator of ios
112  {
113  hasAnyData = false;
114  hasFSData = false;
115  hasXsec = false;
116  // theData.close();
117  return; // no data for exactly this isotope and FS
118  }
119  // here we go
120  G4int infoType, dataType, dummy=INT_MAX;
121  hasFSData = false;
122  while (theData >> infoType)
123  {
124  theData >> dataType;
125 
126  if(dummy==INT_MAX) theData >> dummy >> dummy;
127  if(dataType==3)
128  {
129  G4int total;
130  theData >> total;
131  theXsection->Init(theData, total, CLHEP::eV);
132  }
133  else if(dataType==4)
134  {
136  theAngularDistribution->Init(theData);
137  hasFSData = true;
138  }
139  else if(dataType==5)
140  {
142  theEnergyDistribution->Init(theData);
143  hasFSData = true;
144  }
145  else if(dataType==6)
146  {
148  //- G4cout << this << " BaseFS theEnergyAngData " << theEnergyAngData << G4endl; //GDEB
149  theEnergyAngData->Init(theData);
150  hasFSData = true;
151  }
152  else if(dataType==12)
153  {
155  theFinalStatePhotons->InitMean(theData);
156  hasFSData = true;
157  }
158  else if(dataType==13)
159  {
162  hasFSData = true;
163  }
164  else if(dataType==14)
165  {
167  hasFSData = true;
168  }
169  else if(dataType==15)
170  {
172  hasFSData = true;
173  }
174  else
175  {
176  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4ParticleHPInelasticBaseFS");
177  }
178  }
179  // theData.close();
180 }
181 
183  G4ParticleDefinition ** theDefs,
184  G4int nDef)
185 {
186 
187 // prepare neutron
188  theResult.Clear();
189  G4double eKinetic = theTrack.GetKineticEnergy();
190  const G4HadProjectile *hadProjectile = &theTrack;
191  G4ReactionProduct incidReactionProduct( const_cast<G4ParticleDefinition *>(hadProjectile->GetDefinition()) );
192  incidReactionProduct.SetMomentum( hadProjectile->Get4Momentum().vect() );
193  incidReactionProduct.SetKineticEnergy( eKinetic );
194 
195 // prepare target
196  G4double targetMass;
197  G4double eps = 0.0001;
198  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
200 
201  if(theEnergyAngData!=0)
202  { targetMass = theEnergyAngData->GetTargetMass(); }
204  { targetMass = theAngularDistribution->GetTargetMass(); }
205 
206 //080731a
207 //110512 TKDB ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded.
208 //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;
209  if ( targetMass == 0 )
210  {
211  //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;
212  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / theProjectile->GetPDGMass();
213  }
214 
215  G4Nucleus aNucleus;
217  G4ThreeVector neuVelo = (1./hadProjectile->GetDefinition()->GetPDGMass())*incidReactionProduct.GetMomentum();
218  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
219 
220  theTarget.SetDefinition( G4IonTable::GetIonTable()->GetIon( G4int(theBaseZ), G4int(theBaseA) , 0.0 ) );
221 
222 // prepare energy in target rest frame
223  G4ReactionProduct boosted;
224  boosted.Lorentz(incidReactionProduct, theTarget);
225  eKinetic = boosted.GetKineticEnergy();
226  G4double orgMomentum = boosted.GetMomentum().mag();
227 
228 // Take N-body phase-space distribution, if no other data present.
229  if(!HasFSData()) // adding the residual is trivial here @@@
230  {
231  G4ParticleHPNBodyPhaseSpace thePhaseSpaceDistribution;
232  G4double aPhaseMass=0;
233  G4int ii;
234  for(ii=0; ii<nDef; ii++)
235  {
236  aPhaseMass+=theDefs[ii]->GetPDGMass();
237  }
238  thePhaseSpaceDistribution.Init(aPhaseMass, nDef);
239  thePhaseSpaceDistribution.SetProjectileRP(&incidReactionProduct);
240  thePhaseSpaceDistribution.SetTarget(&theTarget);
241  for(ii=0; ii<nDef; ii++)
242  {
243  G4double massCode = 1000.*std::abs(theDefs[ii]->GetPDGCharge());
244  massCode += theDefs[ii]->GetBaryonNumber();
245  G4double dummy = 0;
246  G4ReactionProduct * aSec = thePhaseSpaceDistribution.Sample(eKinetic, massCode, dummy);
247  aSec->Lorentz(*aSec, -1.*theTarget);
248  G4DynamicParticle * aPart = new G4DynamicParticle();
249  aPart->SetDefinition(aSec->GetDefinition());
250  aPart->SetMomentum(aSec->GetMomentum());
251  delete aSec;
252  theResult.AddSecondary(aPart);
253 #ifdef G4PHPDEBUG
254  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply NoFSData add secondary " << aPart->GetParticleDefinition()->GetParticleName() << " E= " << aPart->GetKineticEnergy() << " NSECO " << theResult.GetNumberOfSecondaries() << G4endl;
255 #endif
256  }
258  //TK120607
259  //Final momentum check should be done before return
261  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
262  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
263  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
264  adjust_final_state ( init_4p_lab );
265 
266  return;
267  }
268 
269 // set target and neutron in the relevant exit channel
270  if(theAngularDistribution!=0)
271  {
272  theAngularDistribution->SetTarget(theTarget);
273  theAngularDistribution->SetProjectileRP(incidReactionProduct);
274  }
275  else if(theEnergyAngData!=0)
276  {
277  theEnergyAngData->SetTarget(theTarget);
278  theEnergyAngData->SetProjectileRP(incidReactionProduct);
279  }
280 
281  G4ReactionProductVector * tmpHadrons = 0;
282  G4int ii = 0, dummy;
283  unsigned int i;
284 
285  if(theEnergyAngData != 0)
286  {
287  tmpHadrons = theEnergyAngData->Sample(eKinetic);
288  }
289  else if(theAngularDistribution!= 0)
290  {
291  G4bool * Done = new G4bool[nDef];
292  G4int i0;
293  for(i0=0; i0<nDef; i0++) Done[i0] = false;
294  if(tmpHadrons == 0)
295  {
296  tmpHadrons = new G4ReactionProductVector;
297  }
298  else
299  {
300  for(i=0; i<tmpHadrons->size(); i++)
301  {
302  for(ii=0; ii<nDef; ii++)
303  if(!Done[ii] && tmpHadrons->operator[](i)->GetDefinition() == theDefs[ii])
304  Done[ii] = true;
305  }
306 #ifdef G4PHPDEBUG
307  if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply secondary previously added " << tmpHadrons->operator[](i)->GetDefinition()->GetParticleName() << " E= " << tmpHadrons->operator[](i)->GetKineticEnergy() << G4endl;
308 #endif
309  }
310  G4ReactionProduct * aHadron;
311  G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)));
312  G4ThreeVector bufferedDirection(0,0,0);
313  for(i0=0; i0<nDef; i0++)
314  {
315  if(!Done[i0])
316  {
317  aHadron = new G4ReactionProduct;
318  if(theEnergyDistribution!=0)
319  {
320  aHadron->SetDefinition(theDefs[i0]);
321  aHadron->SetKineticEnergy(theEnergyDistribution->Sample(eKinetic, dummy));
322  }
323  else if(nDef == 1)
324  {
325  aHadron->SetDefinition(theDefs[i0]);
326  aHadron->SetKineticEnergy(eKinetic);
327  }
328  else if(nDef == 2)
329  {
330  aHadron->SetDefinition(theDefs[i0]);
331  aHadron->SetKineticEnergy(50*CLHEP::MeV);
332  }
333  else
334  {
335  throw G4HadronicException(__FILE__, __LINE__, "No energy distribution to sample from in InelasticBaseFS::BaseApply");
336  }
338  if(theEnergyDistribution==0 && nDef == 2)
339  {
340  if(i0==0)
341  {
342  G4double mass1 = theDefs[0]->GetPDGMass();
343  G4double mass2 = theDefs[1]->GetPDGMass();
344  G4double massn = theProjectile->GetPDGMass();
345  G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge());
346  G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber();
347  G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1);
348  G4double availableEnergy = eKinetic+massn+localMass-mass1-mass2-concreteMass;
349  // available kinetic energy in CMS (non relativistic)
350  G4double emin = availableEnergy+mass1+mass2 - std::sqrt((mass1+mass2)*(mass1+mass2)+orgMomentum*orgMomentum);
351  G4double p1=std::sqrt(2.*mass2*emin);
352  bufferedDirection = p1*aHadron->GetMomentum().unit();
353 #ifdef G4PHPDEBUG
354  if(getenv("G4ParticleHPDebug")) // @@@@@ verify the nucleon counting...
355  {
356  G4cout << "G4ParticleHPInelasticBaseFS "<<z1<<" "<<theBaseZ<<" "<<a1<<" "<<theBaseA<<" "<<availableEnergy<<" "
357  << emin<<G4endl;
358  }
359 #endif
360  }
361  else
362  {
363  bufferedDirection = -bufferedDirection;
364  }
365  // boost from cms to lab
366 #ifdef G4PHPDEBUG
367  if(getenv("G4ParticleHPDebug"))
368  {
369  G4cout << " G4ParticleHPInelasticBaseFS "<<bufferedDirection.mag2()<<G4endl;
370  }
371 #endif
372  aHadron->SetTotalEnergy( std::sqrt(aHadron->GetMass()*aHadron->GetMass()
373  +bufferedDirection.mag2()) );
374  aHadron->SetMomentum(bufferedDirection);
375  aHadron->Lorentz(*aHadron, -1.*(theTarget+incidReactionProduct));
376 #ifdef G4PHPDEBUG
377  if(getenv("G4ParticleHPDebug"))
378  {
379  G4cout << " G4ParticleHPInelasticBaseFS "<<aHadron->GetTotalEnergy()<<" "<<aHadron->GetMomentum()<<G4endl;
380  }
381 #endif
382  }
383  tmpHadrons->push_back(aHadron);
384 #ifdef G4PHPDEBUG
385  if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply FSData add secondary " << aHadron->GetDefinition()->GetParticleName() << " E= " << aHadron->GetKineticEnergy() << G4endl;
386 #endif
387  }
388  }
389  delete [] Done;
390  }
391  else
392  {
393  throw G4HadronicException(__FILE__, __LINE__, "No data to create the neutrons in NInelasticFS");
394  }
395 
396  G4ReactionProductVector * thePhotons = 0;
397  if(theFinalStatePhotons!=0)
398  {
399  // the photon distributions are in the Nucleus rest frame.
400  G4ReactionProduct boosted_tmp;
401  boosted_tmp.Lorentz(incidReactionProduct, theTarget);
402  G4double anEnergy = boosted_tmp.GetKineticEnergy();
403  thePhotons = theFinalStatePhotons->GetPhotons(anEnergy);
404  if(thePhotons!=0)
405  {
406  for(i=0; i<thePhotons->size(); i++)
407  {
408  // back to lab
409  thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), -1.*theTarget);
410  }
411  }
412  }
413  else if(theEnergyAngData!=0)
414  {
415 
416  // PA130927: do not create photons to adjust binding energy
417  G4bool bAdjustPhotons = true;
418 #ifdef PHP_AS_HP
419  bAdjustPhotons = true;
420 #else
421  if ( getenv( "G4PHP_DO_NOT_ADJUST_FINAL_STATE" ) ) bAdjustPhotons = false;
422 #endif
423 
424  if( bAdjustPhotons ) {
425  G4double theGammaEnergy = theEnergyAngData->GetTotalMeanEnergy();
426  G4double anEnergy = boosted.GetKineticEnergy();
427  theGammaEnergy = anEnergy-theGammaEnergy;
428  theGammaEnergy += theNuclearMassDifference;
429  G4double eBindProducts = 0;
430  G4double eBindN = 0;
431  G4double eBindP = 0;
436  G4int ia=0;
437  for(i=0; i<tmpHadrons->size(); i++)
438  {
439  if(tmpHadrons->operator[](i)->GetDefinition() == G4Neutron::Neutron())
440  {
441  eBindProducts+=eBindN;
442  }
443  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Proton::Proton())
444  {
445  eBindProducts+=eBindP;
446  }
447  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Deuteron::Deuteron())
448  {
449  eBindProducts+=eBindD;
450  }
451  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Triton::Triton())
452  {
453  eBindProducts+=eBindT;
454  }
455  else if(tmpHadrons->operator[](i)->GetDefinition() == G4He3::He3())
456  {
457  eBindProducts+=eBindHe3;
458  }
459  else if(tmpHadrons->operator[](i)->GetDefinition() == G4Alpha::Alpha())
460  {
461  eBindProducts+=eBindA;
462  ia++;
463  }
464  }
465 
466  theGammaEnergy += eBindProducts;
467 
468 #ifdef G4PHPDEBUG
469  if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply gamma Energy " << theGammaEnergy << " eBindProducts " << eBindProducts << G4endl;
470 #endif
471 
472  //101111
473  //Special treatment for Be9 + n -> 2n + Be8 -> 2n + a + a
474  if ( (G4int)(theBaseZ+eps) == 4 && (G4int)(theBaseA+eps) == 9 )
475  {
476  // This only valid for G4NDL3.13,,,
477  if ( std::abs( theNuclearMassDifference -
480  && ia == 2 )
481  {
482  theGammaEnergy -= (2*eBindA);
483  }
484  }
485 
486  G4ReactionProductVector * theOtherPhotons = 0;
487  G4int iLevel;
488  while(theGammaEnergy>=theGammas.GetLevelEnergy(0))
489  {
490  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
491  {
492  if(theGammas.GetLevelEnergy(iLevel)<theGammaEnergy) break;
493  }
494  if(iLevel==0||iLevel==theGammas.GetNumberOfLevels()-1)
495  {
496  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
497 #ifdef G4PHPDEBUG
498  if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply adding gamma from level " << iLevel << theOtherPhotons->operator[](ii)->GetKineticEnergy() << G4endl;
499 #endif
500  }
501  else
502  {
503  G4double random = G4UniformRand();
504  G4double eLow = theGammas.GetLevelEnergy(iLevel);
505  G4double eHigh = theGammas.GetLevelEnergy(iLevel+1);
506  if(random > (eHigh-eLow)/(theGammaEnergy-eLow)) iLevel++;
507  theOtherPhotons = theGammas.GetDecayGammas(iLevel);
508  }
509  if(thePhotons==0) thePhotons = new G4ReactionProductVector;
510  if(theOtherPhotons != 0)
511  {
512  for(unsigned int iii=0; iii<theOtherPhotons->size(); iii++)
513  {
514  thePhotons->push_back(theOtherPhotons->operator[](iii));
515 #ifdef G4PHPDEBUG
516  if( getenv("G4ParticleHPDebug"))
517  G4cout << iii << " G4ParticleHPInelasticBaseFS::BaseApply adding gamma " << theOtherPhotons->operator[](iii)->GetKineticEnergy() << G4endl;
518 #endif
519  }
520  delete theOtherPhotons;
521  }
522  theGammaEnergy -= theGammas.GetLevelEnergy(iLevel);
523  if(iLevel == -1) break;
524  }
525  }
526  }
527 
528  // fill the result
529  unsigned int nSecondaries = tmpHadrons->size();
530  unsigned int nPhotons = 0;
531  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
532  nSecondaries += nPhotons;
533  G4DynamicParticle * theSec;
534 #ifdef G4PHPDEBUG
535  if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply N hadrons " << nSecondaries-nPhotons << G4endl;
536 #endif
537 
538  for(i=0; i<nSecondaries-nPhotons; i++)
539  {
540  theSec = new G4DynamicParticle;
541  theSec->SetDefinition(tmpHadrons->operator[](i)->GetDefinition());
542  theSec->SetMomentum(tmpHadrons->operator[](i)->GetMomentum());
543  theResult.AddSecondary(theSec);
544 #ifdef G4PHPDEBUG
545  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply add secondary2 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.GetNumberOfSecondaries() << G4endl;
546 #endif
547  if( getenv("G4PHPTEST") ) G4cout << " InelasticBaseFS COS THETA " << std::cos(theSec->GetMomentum().theta()) << " " << (theSec->GetMomentum().theta()) << " " << theSec->GetMomentum() << " E "<< theSec->GetKineticEnergy() << " " << theSec->GetDefinition()->GetParticleName() << G4endl; //GDEB
548  delete tmpHadrons->operator[](i);
549  }
550 #ifdef G4PHPDEBUG
551  if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply N photons " << nPhotons << G4endl;
552 #endif
553  if(thePhotons != 0)
554  {
555  for(i=0; i<nPhotons; i++)
556  {
557  theSec = new G4DynamicParticle;
558  theSec->SetDefinition(thePhotons->operator[](i)->GetDefinition());
559  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
560  theResult.AddSecondary(theSec);
561 #ifdef G4PHPDEBUG
562  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply add secondary3 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.GetNumberOfSecondaries() << G4endl;
563 #endif
564  delete thePhotons->operator[](i);
565  }
566  }
567 
568 // some garbage collection
569  delete thePhotons;
570  delete tmpHadrons;
571 
572 //080721
574  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
575  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
576  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
577  adjust_final_state ( init_4p_lab );
578 
579 // clean up the primary neutron
581 }
static G4ParticleHPManager * GetInstance()
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass)
G4ReactionProductVector * GetDecayGammas(G4int aLevel)
static const double MeV
Definition: G4SIunits.hh:193
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetMomentum(const G4ThreeVector &momentum)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void SetKineticEnergy(const G4double en)
static const G4double a1
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetProjectileRP(G4ReactionProduct &aIncidentPart)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:463
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
void GetDataStream(G4String, std::istringstream &iss)
G4ParticleHPEnAngCorrelation * theEnergyAngData
static const G4double eps
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double Sample(G4double anEnergy, G4int &it)
G4ReactionProductVector * Sample(G4double anEnergy)
const G4String & GetParticleName() const
void InitPartials(std::istream &aDataFile)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetProjectileRP(const G4ReactionProduct &anIncidentParticleRP)
G4ParticleDefinition * theProjectile
void Init(G4double aMass, G4int aCount)
void BaseApply(const G4HadProjectile &theTrack, G4ParticleDefinition **theDefs, G4int nDef)
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetTarget(const G4ReactionProduct &aTarget)
void adjust_final_state(G4LorentzVector)
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
const G4ParticleDefinition * GetDefinition() const
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
void SampleAndUpdate(G4ReactionProduct &anIncidentParticle)
const G4ParticleDefinition * GetDefinition() const
G4bool InitMean(std::istream &aDataFile)
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
void SetTotalEnergy(const G4double en)
const G4double p1
G4ErrorTarget * theTarget
Definition: errprop.cc:59
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void Init(std::istream &aDataFile)
const G4ParticleDefinition * GetParticleDefinition() const
static const G4double A[nN]
const G4LorentzVector & Get4Momentum() const
void SetProjectileRP(G4ReactionProduct *aIncidentParticleRP)
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
#define INT_MAX
Definition: templates.hh:111
G4double GetKineticEnergy() const
static const double eV
Definition: G4SIunits.hh:194
G4double GetTotalEnergy() const
void SetTarget(G4ReactionProduct *aTarget)
G4double total(Particle const *const p1, Particle const *const p2)
void InitGammas(G4double AR, G4double ZR)
G4double GetPDGMass() const
static G4double GetBindingEnergy(const G4int A, const G4int Z)
G4ParticleHPEnergyDistribution * theEnergyDistribution
void SetTarget(G4ReactionProduct &aTarget)
G4ReactionProductVector * GetPhotons(G4double anEnergy)
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
G4ParticleHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
const G4Material * GetMaterial() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
static const double keV
Definition: G4SIunits.hh:195
void InitEnergies(std::istream &aDataFile)
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &bit, G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ParticleHPPhotonDist * theFinalStatePhotons
void Init(std::istream &aDataFile)
void InitAngular(std::istream &aDataFile)
G4double GetPDGCharge() const
static G4He3 * He3()
Definition: G4He3.cc:94
G4double GetLevelEnergy(G4int aLevel)
G4double GetMass() const
G4int GetNumberOfSecondaries() const
G4ParticleHPAngular * theAngularDistribution
G4ThreeVector GetMomentum() const
void Init(std::istream &aDataFile)
CLHEP::HepLorentzVector G4LorentzVector