Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParticleHPInelasticCompFS Class Referenceabstract

#include <G4ParticleHPInelasticCompFS.hh>

Inheritance diagram for G4ParticleHPInelasticCompFS:
Collaboration diagram for G4ParticleHPInelasticCompFS:

Public Member Functions

 G4ParticleHPInelasticCompFS ()
 
virtual ~G4ParticleHPInelasticCompFS ()
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &aSFType, G4ParticleDefinition *)
 
void InitGammas (G4double AR, G4double ZR)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack)=0
 
virtual G4ParticleHPFinalStateNew ()=0
 
virtual G4double GetXsec (G4double anEnergy)
 
virtual G4ParticleHPVectorGetXsec ()
 
G4int SelectExitChannel (G4double eKinetic)
 
void CompositeApply (const G4HadProjectile &theTrack, G4ParticleDefinition *aHadron)
 
void InitDistributionInitialState (G4ReactionProduct &anIncidentPart, G4ReactionProduct &aTarget, G4int it)
 
- Public Member Functions inherited from G4ParticleHPFinalState
 G4ParticleHPFinalState ()
 
virtual ~G4ParticleHPFinalState ()
 
void Init (G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile)
 
G4bool HasXsec ()
 
G4bool HasFSData ()
 
G4bool HasAnyData ()
 
void SetA_Z (G4double anA, G4double aZ, G4int aM=0)
 
G4double GetZ ()
 
G4double GetN ()
 
G4double GetA ()
 
G4int GetM ()
 
void SetAZMs (G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
 
void SetProjectile (G4ParticleDefinition *projectile)
 

Protected Attributes

G4ParticleHPVectortheXsection [51]
 
G4ParticleHPEnergyDistributiontheEnergyDistribution [51]
 
G4ParticleHPAngulartheAngularDistribution [51]
 
G4ParticleHPEnAngCorrelationtheEnergyAngData [51]
 
G4ParticleHPPhotonDisttheFinalStatePhotons [51]
 
G4ParticleHPDeExGammas theGammas
 
G4String gammaPath
 
std::vector< G4doubleQI
 
std::vector< G4intLR
 
- Protected Attributes inherited from G4ParticleHPFinalState
G4bool hasXsec
 
G4bool hasFSData
 
G4bool hasAnyData
 
G4ParticleHPNames theNames
 
G4Cache< G4HadFinalState * > theResult
 
G4double theBaseA
 
G4double theBaseZ
 
G4int theBaseM
 
G4int theNDLDataZ
 
G4int theNDLDataA
 
G4int theNDLDataM
 
G4ParticleDefinitiontheProjectile
 

Additional Inherited Members

- Protected Member Functions inherited from G4ParticleHPFinalState
void adjust_final_state (G4LorentzVector)
 
G4bool DoNotAdjustFinalState ()
 

Detailed Description

Definition at line 43 of file G4ParticleHPInelasticCompFS.hh.

Constructor & Destructor Documentation

G4ParticleHPInelasticCompFS::G4ParticleHPInelasticCompFS ( )
inline

Definition at line 47 of file G4ParticleHPInelasticCompFS.hh.

48  {
49 
50  QI.resize(51);
51  LR.resize(51);
52  for(G4int i=0; i<51; i++)
53  {
54  hasXsec = true;
55  theXsection[i] = 0;
56  theEnergyDistribution[i] = 0;
58  theEnergyAngData[i] = 0;
59  theFinalStatePhotons[i] = 0;
60  QI[i]=0.0;
61  LR[i]=0;
62  }
63 
64  }
G4ParticleHPEnergyDistribution * theEnergyDistribution[51]
G4ParticleHPPhotonDist * theFinalStatePhotons[51]
int G4int
Definition: G4Types.hh:78
G4ParticleHPEnAngCorrelation * theEnergyAngData[51]
G4ParticleHPAngular * theAngularDistribution[51]
virtual G4ParticleHPInelasticCompFS::~G4ParticleHPInelasticCompFS ( )
inlinevirtual

Definition at line 65 of file G4ParticleHPInelasticCompFS.hh.

66  {
67  for(G4int i=0; i<51; i++)
68  {
69  if(theXsection[i] != 0) delete theXsection[i];
70  if(theEnergyDistribution[i] != 0) delete theEnergyDistribution[i];
71  if(theAngularDistribution[i] != 0) delete theAngularDistribution[i];
72  if(theEnergyAngData[i] != 0) delete theEnergyAngData[i];
73  if(theFinalStatePhotons[i] != 0) delete theFinalStatePhotons[i];
74  }
75  }
G4ParticleHPEnergyDistribution * theEnergyDistribution[51]
G4ParticleHPPhotonDist * theFinalStatePhotons[51]
int G4int
Definition: G4Types.hh:78
G4ParticleHPEnAngCorrelation * theEnergyAngData[51]
G4ParticleHPAngular * theAngularDistribution[51]

Member Function Documentation

virtual G4HadFinalState* G4ParticleHPInelasticCompFS::ApplyYourself ( const G4HadProjectile theTrack)
pure virtual
void G4ParticleHPInelasticCompFS::CompositeApply ( const G4HadProjectile theTrack,
G4ParticleDefinition aHadron 
)

Definition at line 230 of file G4ParticleHPInelasticCompFS.cc.

231 {
232 
233 // prepare neutron
234  if ( theResult.Get() == NULL ) theResult.Put( new G4HadFinalState );
235  theResult.Get()->Clear();
236  G4double eKinetic = theTrack.GetKineticEnergy();
237  const G4HadProjectile *hadProjectile = &theTrack;
238  G4ReactionProduct incidReactionProduct( const_cast<G4ParticleDefinition *>(hadProjectile->GetDefinition()) ); // incidReactionProduct
239  incidReactionProduct.SetMomentum( hadProjectile->Get4Momentum().vect() );
240  incidReactionProduct.SetKineticEnergy( eKinetic );
241 
242 // prepare target
243  G4int i;
244  for(i=0; i<50; i++)
245  { if(theXsection[i] != 0) { break; } }
246 
247  G4double targetMass=0;
248  G4double eps = 0.0001;
249  targetMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
250 #ifdef G4PHPDEBUG
251  if( getenv("G4ParticleHPDebug")) G4cout <<this <<" G4ParticleHPInelasticCompFS::CompositeApply A " <<theBaseA <<" Z " <<theBaseZ <<" incident " <<hadProjectile->GetDefinition()->GetParticleName() <<G4endl;
252 #endif
253 // if(theEnergyAngData[i]!=0)
254 // targetMass = theEnergyAngData[i]->GetTargetMass();
255 // else if(theAngularDistribution[i]!=0)
256 // targetMass = theAngularDistribution[i]->GetTargetMass();
257 // else if(theFinalStatePhotons[50]!=0)
258 // targetMass = theFinalStatePhotons[50]->GetTargetMass();
260  G4Nucleus aNucleus;
261  //G4ThreeVector neuVelo = (1./hadProjectile->GetDefinition()->GetPDGMass())*incidReactionProduct.GetMomentum();
262  //theTarget = aNucleus.GetBiasedThermalNucleus( targetMass/hadProjectile->GetDefinition()->GetPDGMass() , neuVelo, theTrack.GetMaterial()->GetTemperature());
263  //G4Nucleus::GetBiasedThermalNucleus requests normalization of mass and velocity in neutron mass
264  G4ThreeVector neuVelo = ( 1./G4Neutron::Neutron()->GetPDGMass() )*incidReactionProduct.GetMomentum();
265  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass/G4Neutron::Neutron()->GetPDGMass()
266  , neuVelo, theTrack.GetMaterial()->GetTemperature() );
267 
268  theTarget.SetDefinition( G4IonTable::GetIonTable()->GetIon( G4int(theBaseZ), G4int(theBaseA) , 0.0 ) ); //XX
269 
270 // prepare the residual mass
271  G4double residualMass=0;
272  G4double residualZ = theBaseZ + theProjectile->GetPDGCharge() - aDefinition->GetPDGCharge();
273  G4double residualA = theBaseA + theProjectile->GetBaryonNumber() - aDefinition->GetBaryonNumber();
274  residualMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(residualA+eps), static_cast<G4int>(residualZ+eps));
275 
276 // prepare energy in target rest frame
277  G4ReactionProduct boosted;
278  boosted.Lorentz(incidReactionProduct, theTarget);
279  eKinetic = boosted.GetKineticEnergy();
280 // G4double momentumInCMS = boosted.GetTotalMomentum();
281 
282 // select exit channel for composite FS class.
283  G4int it = SelectExitChannel( eKinetic );
284 
285 // set target and neutron in the relevant exit channel
286  InitDistributionInitialState(incidReactionProduct, theTarget, it);
287 
288  G4ReactionProductVector * thePhotons = 0;
289  G4ReactionProductVector * theParticles = 0;
290  G4ReactionProduct aHadron;
291  aHadron.SetDefinition(aDefinition); // what if only cross-sections exist ==> Na 23 11 @@@@
292  G4double availableEnergy = incidReactionProduct.GetKineticEnergy() + incidReactionProduct.GetMass() - aHadron.GetMass() +
293  (targetMass - residualMass);
294 //080730c
295  if ( availableEnergy < 0 )
296  {
297  //G4cout << "080730c Adjust availavleEnergy " << G4endl;
298  availableEnergy = 0;
299  }
300  G4int nothingWasKnownOnHadron = 0;
301  G4int dummy;
302  G4double eGamm = 0;
303  G4int iLevel=it-1;
304 
305 // TK without photon has it = 0
306  if( 50 == it )
307  {
308 
309 // TK Excitation level is not determined
310  iLevel=-1;
311  aHadron.SetKineticEnergy(availableEnergy*residualMass/
312  (aHadron.GetMass()+residualMass));
313 
314  //aHadron.SetMomentum(incidReactionProduct.GetMomentum()*(1./incidReactionProduct.GetTotalMomentum())*
315  // std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
316  // aHadron.GetMass()*aHadron.GetMass()));
317 
318  //TK add safty 100909
319  G4double p2 = ( aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy() - aHadron.GetMass()*aHadron.GetMass() );
320  G4double p = 0.0;
321  if ( p2 > 0.0 ) p = std::sqrt( p );
322 
323  aHadron.SetMomentum(incidReactionProduct.GetMomentum()*(1./incidReactionProduct.GetTotalMomentum())*p );
324 
325  }
326  else
327  {
328  while ( iLevel!=-1 && theGammas.GetLevel(iLevel) == 0 ) { iLevel--; } // Loop checking, 11.05.2015, T. Koi
329  }
330 
331 
332  if ( theAngularDistribution[it] != 0 ) // MF4
333  {
334  if(theEnergyDistribution[it]!=0) // MF5
335  {
336  //************************************************************
337  /*
338  aHadron.SetKineticEnergy(theEnergyDistribution[it]->Sample(eKinetic, dummy));
339  G4double eSecN = aHadron.GetKineticEnergy();
340  */
341  //************************************************************
342  //EMendoza --> maximum allowable energy should be taken into account.
343  G4double dqi = 0.0;
344  if ( QI[it] < 0 || 849 < QI[it] ) dqi = QI[it]; //For backword compatibility QI introduced since G4NDL3.15
345  G4double MaxEne=eKinetic+dqi;
346  G4double eSecN;
347 
348  G4int icounter=0;
349  G4int icounter_max=1024;
350  do {
351  icounter++;
352  if ( icounter > icounter_max ) {
353  G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
354  break;
355  }
356  eSecN=theEnergyDistribution[it]->Sample(eKinetic, dummy);
357  }while(eSecN>MaxEne); // Loop checking, 11.05.2015, T. Koi
358  aHadron.SetKineticEnergy(eSecN);
359  //************************************************************
360  eGamm = eKinetic-eSecN;
361  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
362  {
363  if(theGammas.GetLevelEnergy(iLevel)<eGamm) break;
364  }
365  G4double random = 2*G4UniformRand();
366  iLevel+=G4int(random);
367  if(iLevel>theGammas.GetNumberOfLevels()-1)iLevel = theGammas.GetNumberOfLevels()-1;
368  }
369  else
370  {
371  G4double eExcitation = 0;
372  if(iLevel>=0) eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();
373  while (eKinetic-eExcitation < 0 && iLevel>0) // Loop checking, 11.05.2015, T. Koi
374  {
375  iLevel--;
376  eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();
377  }
378  //110610TK BEGIN
379  //Use QI value for calculating excitation energy of residual.
380  G4bool useQI=false;
381  G4double dqi = QI[it];
382  if ( dqi < 0 || 849 < dqi ) useQI = true; //Former libraies does not have values of this range
383 
384  if ( useQI )
385  {
386  // QI introudced since G4NDL3.15
387  eExcitation = -QI[it];
388  //Re-evluate iLevel based on this eExcitation
389  iLevel = 0;
390  G4bool find = false;
391  G4int imaxEx = 0;
392  while( theGammas.GetLevel(iLevel+1) != 0 ) // Loop checking, 11.05.2015, T. Koi
393  {
394  G4double maxEx = 0.0;
395  if ( maxEx < theGammas.GetLevel(iLevel)->GetLevelEnergy() )
396  {
397  maxEx = theGammas.GetLevel(iLevel)->GetLevelEnergy();
398  imaxEx = iLevel;
399  }
400  if ( eExcitation < theGammas.GetLevel(iLevel)->GetLevelEnergy() )
401  {
402  find = true;
403  iLevel--;
404  // very small eExcitation, iLevel becomes -1, this is protected below.
405  if ( iLevel == -1 ) iLevel = 0; // But cause energy trouble.
406  break;
407  }
408  iLevel++;
409  }
410  // In case, cannot find proper level, then use the maximum level.
411  if ( !find ) iLevel = imaxEx;
412  }
413  //110610TK END
414 
415  if(getenv("G4ParticleHPDebug") && eKinetic-eExcitation < 0)
416  {
417  throw G4HadronicException(__FILE__, __LINE__, "SEVERE: InelasticCompFS: Consistency of data not good enough, please file report");
418  }
419  if(eKinetic-eExcitation < 0) eExcitation = 0;
420  if(iLevel!= -1) aHadron.SetKineticEnergy(eKinetic - eExcitation);
421 
422  }
424 
425  if( theFinalStatePhotons[it] == 0 )
426  {
427  //G4cout << "110610 USE Gamma Level" << G4endl;
428 // TK comment Most n,n* eneter to this
429  thePhotons = theGammas.GetDecayGammas(iLevel);
430  eGamm -= theGammas.GetLevelEnergy(iLevel);
431  if(eGamm>0) // @ ok for now, but really needs an efficient way of correllated sampling @
432  {
433  G4ReactionProduct * theRestEnergy = new G4ReactionProduct;
434  theRestEnergy->SetDefinition(G4Gamma::Gamma());
435  theRestEnergy->SetKineticEnergy(eGamm);
436  G4double costh = 2.*G4UniformRand()-1.;
438  theRestEnergy->SetMomentum(eGamm*std::sin(std::acos(costh))*std::cos(phi),
439  eGamm*std::sin(std::acos(costh))*std::sin(phi),
440  eGamm*costh);
441  if(thePhotons == 0) { thePhotons = new G4ReactionProductVector; }
442  thePhotons->push_back(theRestEnergy);
443  }
444  }
445  }
446  else if(theEnergyAngData[it] != 0) // MF6
447  {
448 
449  theParticles = theEnergyAngData[it]->Sample(eKinetic);
450 
451  //141017 Fix BEGIN
452  //Adjust A and Z in the case of miss much between selected data and target nucleus
453  if ( theParticles != NULL ) {
454  G4int sumA = 0;
455  G4int sumZ = 0;
456  G4int maxA = 0;
457  G4int jAtMaxA = 0;
458  for ( G4int j = 0 ; j != (G4int)theParticles->size() ; j++ ) {
459  if ( theParticles->at(j)->GetDefinition()->GetBaryonNumber() > maxA ) {
460  maxA = theParticles->at(j)->GetDefinition()->GetBaryonNumber();
461  jAtMaxA = j;
462  }
463  sumA += theParticles->at(j)->GetDefinition()->GetBaryonNumber();
464  sumZ += G4int( theParticles->at(j)->GetDefinition()->GetPDGCharge() + eps );
465  }
466  G4int dA = (G4int)theBaseA + hadProjectile->GetDefinition()->GetBaryonNumber() - sumA;
467  G4int dZ = (G4int)theBaseZ + G4int( hadProjectile->GetDefinition()->GetPDGCharge() + eps ) - sumZ;
468  if ( dA < 0 || dZ < 0 ) {
469  G4int newA = theParticles->at(jAtMaxA)->GetDefinition()->GetBaryonNumber() + dA ;
470  G4int newZ = G4int( theParticles->at(jAtMaxA)->GetDefinition()->GetPDGCharge() + eps ) + dZ;
471  G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon ( newZ , newA );
472  theParticles->at( jAtMaxA )->SetDefinition( pd );
473  }
474  }
475  //141017 Fix END
476 
477  }
478  else
479  {
480  // @@@ what to do, if we have photon data, but no info on the hadron itself
481  nothingWasKnownOnHadron = 1;
482  }
483 
484  //G4cout << "theFinalStatePhotons it " << it << G4endl;
485  //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
486  //G4cout << "theFinalStatePhotons it " << it << G4endl;
487  //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
488  //G4cout << "thePhotons " << thePhotons << G4endl;
489 
490  if ( theFinalStatePhotons[it] != 0 )
491  {
492  // the photon distributions are in the Nucleus rest frame.
493  // TK residual rest frame
494  G4ReactionProduct boosted_tmp;
495  boosted_tmp.Lorentz(incidReactionProduct, theTarget);
496  G4double anEnergy = boosted_tmp.GetKineticEnergy();
497  thePhotons = theFinalStatePhotons[it]->GetPhotons(anEnergy);
498  G4double aBaseEnergy = theFinalStatePhotons[it]->GetLevelEnergy();
499  G4double testEnergy = 0;
500  if(thePhotons!=0 && thePhotons->size()!=0)
501  { aBaseEnergy-=thePhotons->operator[](0)->GetTotalEnergy(); }
502  if(theFinalStatePhotons[it]->NeedsCascade())
503  {
504  while(aBaseEnergy>0.01*CLHEP::keV) // Loop checking, 11.05.2015, T. Koi
505  {
506  // cascade down the levels
507  G4bool foundMatchingLevel = false;
508  G4int closest = 2;
509  G4double deltaEold = -1;
510  for(G4int j=1; j<it; j++)
511  {
512  if(theFinalStatePhotons[j]!=0)
513  {
514  testEnergy = theFinalStatePhotons[j]->GetLevelEnergy();
515  }
516  else
517  {
518  testEnergy = 0;
519  }
520  G4double deltaE = std::abs(testEnergy-aBaseEnergy);
521  if(deltaE<0.1*CLHEP::keV)
522  {
523  G4ReactionProductVector * theNext =
524  theFinalStatePhotons[j]->GetPhotons(anEnergy);
525  if ( thePhotons != NULL ) thePhotons->push_back(theNext->operator[](0));
526  aBaseEnergy = testEnergy-theNext->operator[](0)->GetTotalEnergy();
527  delete theNext;
528  foundMatchingLevel = true;
529  break; // ===>
530  }
531  if(theFinalStatePhotons[j]!=0 && ( deltaE<deltaEold||deltaEold<0.) )
532  {
533  closest = j;
534  deltaEold = deltaE;
535  }
536  } // <=== the break goes here.
537  if(!foundMatchingLevel)
538  {
539  G4ReactionProductVector * theNext =
540  theFinalStatePhotons[closest]->GetPhotons(anEnergy);
541  if ( thePhotons != NULL ) thePhotons->push_back(theNext->operator[](0));
542  aBaseEnergy = aBaseEnergy-theNext->operator[](0)->GetTotalEnergy();
543  delete theNext;
544  }
545  }
546  }
547  }
548  unsigned int i0;
549  if(thePhotons!=0)
550  {
551  for(i0=0; i0<thePhotons->size(); i0++)
552  {
553  // back to lab
554  thePhotons->operator[](i0)->Lorentz(*(thePhotons->operator[](i0)), -1.*theTarget);
555  }
556  }
557  //G4cout << "nothingWasKnownOnHadron " << nothingWasKnownOnHadron << G4endl;
558  if(nothingWasKnownOnHadron)
559  {
560 // TKDB 100405
561 // In this case, hadron should be isotropic in CM
562 // mu and p should be correlated
563 //
564  G4double totalPhotonEnergy = 0.0;
565  if ( thePhotons != 0 )
566  {
567  unsigned int nPhotons = thePhotons->size();
568  unsigned int ii0;
569  for ( ii0=0; ii0<nPhotons; ii0++)
570  {
571  //thePhotons has energies at LAB system
572  totalPhotonEnergy += thePhotons->operator[](ii0)->GetTotalEnergy();
573  }
574  }
575 
576  //isotropic distribution in CM
577  G4double mu = 1.0 - 2 * G4UniformRand();
578 
579  // need momentums in target rest frame;
580  G4LorentzVector target_in_LAB ( theTarget.GetMomentum() , theTarget.GetTotalEnergy() );
581  G4ThreeVector boostToTargetRest = -target_in_LAB.boostVector();
582  G4LorentzVector proj_in_LAB = hadProjectile->Get4Momentum();
583 
584  G4DynamicParticle* proj = new G4DynamicParticle( theProjectile , proj_in_LAB.boost( boostToTargetRest ) );
585  G4DynamicParticle* targ = new G4DynamicParticle( G4IonTable::GetIonTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , totalPhotonEnergy ) , G4ThreeVector(0) );
586  G4DynamicParticle* hadron = new G4DynamicParticle( aHadron.GetDefinition() , G4ThreeVector(0) ); // will be fill momentum
587 
588  two_body_reaction ( proj , targ , hadron , mu );
589 
590  G4LorentzVector hadron_in_trag_rest = hadron->Get4Momentum();
591  G4LorentzVector hadron_in_LAB = hadron_in_trag_rest.boost ( -boostToTargetRest );
592  aHadron.SetMomentum( hadron_in_LAB.v() );
593  aHadron.SetKineticEnergy ( hadron_in_LAB.e() - hadron_in_LAB.m() );
594 
595  delete proj;
596  delete targ;
597  delete hadron;
598 
599 //TKDB 100405
600 /*
601  G4double totalPhotonEnergy = 0;
602  if(thePhotons!=0)
603  {
604  unsigned int nPhotons = thePhotons->size();
605  unsigned int i0;
606  for(i0=0; i0<nPhotons; i0++)
607  {
608  totalPhotonEnergy += thePhotons->operator[](i0)->GetTotalEnergy();
609  }
610  }
611  availableEnergy -= totalPhotonEnergy;
612  residualMass += totalPhotonEnergy/theProjectile->GetPDGMass();
613  aHadron.SetKineticEnergy(availableEnergy*residualMass*theProjectile->GetPDGMass()/
614  (aHadron.GetMass()+residualMass*theProjectile->GetPDGMass()));
615  G4double CosTheta = 1.0 - 2.0*G4UniformRand();
616  G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
617  G4double Phi = twopi*G4UniformRand();
618  G4ThreeVector Vector(std::cos(Phi)*SinTheta, std::sin(Phi)*SinTheta, CosTheta);
619  //aHadron.SetMomentum(Vector* std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
620  // aHadron.GetMass()*aHadron.GetMass()));
621  G4double p2 = aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()- aHadron.GetMass()*aHadron.GetMass();
622 
623  G4double p = 0.0;
624  if ( p2 > 0.0 )
625  p = std::sqrt ( p2 );
626 
627  aHadron.SetMomentum( Vector*p );
628 */
629 
630  }
631 
632 // fill the result
633 // Beware - the recoil is not necessarily in the particles...
634 // Can be calculated from momentum conservation?
635 // The idea is that the particles ar emitted forst, and the gammas only once the
636 // recoil is on the residual; assumption is that gammas do not contribute to
637 // the recoil.
638 // This needs more design @@@
639 
640  G4int nSecondaries = 2; // the hadron and the recoil
641  G4bool needsSeparateRecoil = false;
642  G4int totalBaryonNumber = 0;
643  G4int totalCharge = 0;
644  G4ThreeVector totalMomentum(0);
645  if(theParticles != 0)
646  {
647  nSecondaries = theParticles->size();
648  const G4ParticleDefinition * aDef;
649  unsigned int ii0;
650  for(ii0=0; ii0<theParticles->size(); ii0++)
651  {
652  aDef = theParticles->operator[](ii0)->GetDefinition();
653  totalBaryonNumber+=aDef->GetBaryonNumber();
654  totalCharge+=G4int(aDef->GetPDGCharge()+eps);
655  totalMomentum += theParticles->operator[](ii0)->GetMomentum();
656  }
657  if(totalBaryonNumber!=G4int(theBaseA+eps+hadProjectile->GetDefinition()->GetBaryonNumber()))
658  {
659  needsSeparateRecoil = true;
660  nSecondaries++;
661  residualA = G4int(theBaseA+eps+hadProjectile->GetDefinition()->GetBaryonNumber()
662  -totalBaryonNumber);
663  residualZ = G4int(theBaseZ+eps+hadProjectile->GetDefinition()->GetPDGCharge()
664  -totalCharge);
665  }
666  }
667 
668  G4int nPhotons = 0;
669  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
670  nSecondaries += nPhotons;
671 
672  G4DynamicParticle * theSec;
673 
674  if( theParticles==0 )
675  {
676  theSec = new G4DynamicParticle;
677  theSec->SetDefinition(aHadron.GetDefinition());
678  theSec->SetMomentum(aHadron.GetMomentum());
679  theResult.Get()->AddSecondary(theSec);
680 #ifdef G4PHPDEBUG
681  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary1 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
682 #endif
683 
684  aHadron.Lorentz(aHadron, theTarget);
685  G4ReactionProduct theResidual;
687  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0));
688  theResidual.SetKineticEnergy(aHadron.GetKineticEnergy()*aHadron.GetMass()/theResidual.GetMass());
689 
690  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #6
691  //theResidual.SetMomentum(-1.*aHadron.GetMomentum());
692  G4ThreeVector incidentNeutronMomentum = incidReactionProduct.GetMomentum();
693  theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
694 
695  theResidual.Lorentz(theResidual, -1.*theTarget);
696  G4ThreeVector totalPhotonMomentum(0,0,0);
697  if(thePhotons!=0)
698  {
699  for(i=0; i<nPhotons; i++)
700  {
701  totalPhotonMomentum += thePhotons->operator[](i)->GetMomentum();
702  }
703  }
704  theSec = new G4DynamicParticle;
705  theSec->SetDefinition(theResidual.GetDefinition());
706  theSec->SetMomentum(theResidual.GetMomentum()-totalPhotonMomentum);
707  theResult.Get()->AddSecondary(theSec);
708 #ifdef G4PHPDEBUG
709  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary2 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
710 #endif
711  }
712  else
713  {
714  for(i0=0; i0<theParticles->size(); i0++)
715  {
716  theSec = new G4DynamicParticle;
717  theSec->SetDefinition(theParticles->operator[](i0)->GetDefinition());
718  theSec->SetMomentum(theParticles->operator[](i0)->GetMomentum());
719  theResult.Get()->AddSecondary(theSec);
720 #ifdef G4PHPDEBUG
721  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary3 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
722 #endif
723  delete theParticles->operator[](i0);
724  }
725  delete theParticles;
726  if(needsSeparateRecoil && residualZ!=0)
727  {
728  G4ReactionProduct theResidual;
730  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0));
731  G4double resiualKineticEnergy = theResidual.GetMass()*theResidual.GetMass();
732  resiualKineticEnergy += totalMomentum*totalMomentum;
733  resiualKineticEnergy = std::sqrt(resiualKineticEnergy) - theResidual.GetMass();
734 // cout << "Kinetic energy of the residual = "<<resiualKineticEnergy<<endl;
735  theResidual.SetKineticEnergy(resiualKineticEnergy);
736 
737  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #4
738  //theResidual.SetMomentum(-1.*totalMomentum);
739  //G4ThreeVector incidentNeutronMomentum = incidReactionProduct.GetMomentum();
740  //theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
741 //080717 TK Comment still do NOT include photon's mometum which produce by thePhotons
742  theResidual.SetMomentum( incidReactionProduct.GetMomentum() + theTarget.GetMomentum() - totalMomentum );
743 
744  theSec = new G4DynamicParticle;
745  theSec->SetDefinition(theResidual.GetDefinition());
746  theSec->SetMomentum(theResidual.GetMomentum());
747  theResult.Get()->AddSecondary(theSec);
748 #ifdef G4PHPDEBUG
749  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary4 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
750 #endif
751 
752  }
753  }
754  if(thePhotons!=0)
755  {
756  for(i=0; i<nPhotons; i++)
757  {
758  theSec = new G4DynamicParticle;
759  //Bug reported Chao Zhang (Chao.Zhang@usd.edu), Dongming Mei(Dongming.Mei@usd.edu) Feb. 25, 2009
760  //theSec->SetDefinition(G4Gamma::Gamma());
761  theSec->SetDefinition( thePhotons->operator[](i)->GetDefinition() );
762  //But never cause real effect at least with G4NDL3.13 TK
763  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
764  theResult.Get()->AddSecondary(theSec);
765 #ifdef G4PHPDEBUG
766  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary5 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
767 #endif
768 
769  delete thePhotons->operator[](i);
770  }
771 // some garbage collection
772  delete thePhotons;
773  }
774 
775 //080721
776  G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , 0.0 );
777  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
778  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
779  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
780  adjust_final_state ( init_4p_lab );
781 
782 // clean up the primary neutron
784 }
G4ReactionProductVector * GetDecayGammas(G4int aLevel)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ParticleHPEnergyDistribution * theEnergyDistribution[51]
G4Cache< G4HadFinalState * > theResult
void SetMomentum(const G4ThreeVector &momentum)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
static constexpr double keV
G4ParticleHPPhotonDist * theFinalStatePhotons[51]
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
value_type & Get() const
Definition: G4Cache.hh:282
Hep3Vector v() const
static const G4double eps
int G4int
Definition: G4Types.hh:78
G4double Sample(G4double anEnergy, G4int &it)
G4ReactionProductVector * Sample(G4double anEnergy)
const G4String & GetParticleName() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ParticleDefinition * theProjectile
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void adjust_final_state(G4LorentzVector)
const G4ParticleDefinition * GetDefinition() const
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:113
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4ParticleHPEnAngCorrelation * theEnergyAngData[51]
void SampleAndUpdate(G4ReactionProduct &anIncidentParticle)
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
HepLorentzVector & boost(double, double, double)
G4double GetKineticEnergy() const
G4ErrorTarget * theTarget
Definition: errprop.cc:59
G4ParticleHPLevel * GetLevel(G4int i)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4ParticleDefinition * GetParticleDefinition() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void InitDistributionInitialState(G4ReactionProduct &anIncidentPart, G4ReactionProduct &aTarget, G4int it)
const G4LorentzVector & Get4Momentum() const
G4LorentzVector Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4double GetPDGMass() const
double mag2() const
G4ReactionProductVector * GetPhotons(G4double anEnergy)
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
G4ParticleHPAngular * theAngularDistribution[51]
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4int SelectExitChannel(G4double eKinetic)
G4double GetPDGCharge() const
void Put(const value_type &val) const
Definition: G4Cache.hh:286
G4double GetLevelEnergy(G4int aLevel)
G4double GetMass() const
static constexpr double twopi
Definition: SystemOfUnits.h:55
G4int GetNumberOfSecondaries() const

Here is the caller graph for this function:

virtual G4double G4ParticleHPInelasticCompFS::GetXsec ( G4double  anEnergy)
inlinevirtual

Reimplemented from G4ParticleHPFinalState.

Definition at line 80 of file G4ParticleHPInelasticCompFS.hh.

81  {
82  return std::max(0., theXsection[50]->GetY(anEnergy));
83  }
T max(const T t1, const T t2)
brief Return the largest of the two arguments

Here is the call graph for this function:

virtual G4ParticleHPVector* G4ParticleHPInelasticCompFS::GetXsec ( )
inlinevirtual

Reimplemented from G4ParticleHPFinalState.

Definition at line 84 of file G4ParticleHPInelasticCompFS.hh.

84 { return theXsection[50]; }

Here is the caller graph for this function:

void G4ParticleHPInelasticCompFS::Init ( G4double  A,
G4double  Z,
G4int  M,
G4String dirName,
G4String aSFType,
G4ParticleDefinition  
)
virtual

Implements G4ParticleHPFinalState.

Reimplemented in G4ParticleHPNInelasticFS, G4ParticleHPPInelasticFS, and G4ParticleHPTInelasticFS.

Definition at line 81 of file G4ParticleHPInelasticCompFS.cc.

82 {
83  gammaPath = "/Inelastic/Gammas/"; //only in neutron data base
84  if(!getenv("G4NEUTRONHPDATA"))
85  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files where Inelastic/Gammas data is found.");
86  G4String tBase = getenv("G4NEUTRONHPDATA");
87  gammaPath = tBase+gammaPath;
88  G4String tString = dirName;
89  G4bool dbool;
90  G4ParticleHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, tString, aFSType, dbool);
91  G4String filename = aFile.GetName();
92 #ifdef G4PHPDEBUG
93  if( getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPInelasticCompFS::Init FILE " << filename << G4endl;
94 #endif
95 
96  SetAZMs( A, Z, M, aFile );
97  //theBaseA = aFile.GetA();
98  //theBaseZ = aFile.GetZ();
99  //theNDLDataA = (int)aFile.GetA();
100  //theNDLDataZ = aFile.GetZ();
101  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
102  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
103  {
104 #ifdef G4PHPDEBUG
105  if(getenv("G4ParticleHPDebug_NamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
106 #endif
107  hasAnyData = false;
108  hasFSData = false;
109  hasXsec = false;
110  return;
111  }
112  // theBaseA = A;
113  // theBaseZ = G4int(Z+.5);
114 //std::ifstream theData(filename, std::ios::in);
115  std::istringstream theData(std::ios::in);
116  G4ParticleHPManager::GetInstance()->GetDataStream(filename,theData);
117  if(!theData) //"!" is a operator of ios
118  {
119  hasAnyData = false;
120  hasFSData = false;
121  hasXsec = false;
122  // theData.close();
123  return;
124  }
125  // here we go
126  G4int infoType, dataType, dummy;
127  G4int sfType, it;
128  hasFSData = false;
129  while (theData >> infoType) // Loop checking, 11.05.2015, T. Koi
130  {
131  hasFSData = true;
132  theData >> dataType;
133  theData >> sfType >> dummy;
134  it = 50;
135  if(sfType>=600||(sfType<100&&sfType>=50)) it = sfType%50;
136  if(dataType==3)
137  {
138  //theData >> dummy >> dummy;
139  //TK110430
140  // QI and LR introudced since G4NDL3.15
141  G4double dqi;
142  G4int ilr;
143  theData >> dqi >> ilr;
144 
145  QI[ it ] = dqi*CLHEP::eV;
146  LR[ it ] = ilr;
148  G4int total;
149  theData >> total;
150  theXsection[it]->Init(theData, total, CLHEP::eV);
151  //std::cout << theXsection[it]->GetXsec(1*MeV) << std::endl;
152  }
153  else if(dataType==4)
154  {
156  theAngularDistribution[it]->Init(theData);
157  }
158  else if(dataType==5)
159  {
161  theEnergyDistribution[it]->Init(theData);
162  }
163  else if(dataType==6)
164  {
166  // G4cout << this << " CompFS theEnergyAngData " << it << theEnergyAngData[it] << G4endl; //GDEB
167  theEnergyAngData[it]->Init(theData);
168  }
169  else if(dataType==12)
170  {
172  theFinalStatePhotons[it]->InitMean(theData);
173  }
174  else if(dataType==13)
175  {
177  theFinalStatePhotons[it]->InitPartials(theData);
178  }
179  else if(dataType==14)
180  {
181  theFinalStatePhotons[it]->InitAngular(theData);
182  }
183  else if(dataType==15)
184  {
185  theFinalStatePhotons[it]->InitEnergies(theData);
186  }
187  else
188  {
189  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4ParticleHPInelasticCompFS");
190  }
191  }
192  // theData.close();
193 }
static G4ParticleHPManager * GetInstance()
G4ParticleHPEnergyDistribution * theEnergyDistribution[51]
G4ParticleHPPhotonDist * theFinalStatePhotons[51]
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
void GetDataStream(G4String, std::istringstream &iss)
int G4int
Definition: G4Types.hh:78
void InitPartials(std::istream &aDataFile)
G4ParticleDefinition * theProjectile
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
G4GLOB_DLL std::ostream G4cout
G4ParticleHPEnAngCorrelation * theEnergyAngData[51]
double A(double temperature)
G4bool InitMean(std::istream &aDataFile)
bool G4bool
Definition: G4Types.hh:79
static constexpr double eV
G4double total(Particle const *const p1, Particle const *const p2)
#define G4endl
Definition: G4ios.hh:61
G4ParticleHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
void InitEnergies(std::istream &aDataFile)
G4ParticleHPAngular * theAngularDistribution[51]
double G4double
Definition: G4Types.hh:76
void Init(std::istream &aDataFile)
void InitAngular(std::istream &aDataFile)
void Init(std::istream &aDataFile)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleHPInelasticCompFS::InitDistributionInitialState ( G4ReactionProduct anIncidentPart,
G4ReactionProduct aTarget,
G4int  it 
)
inline

Definition at line 87 of file G4ParticleHPInelasticCompFS.hh.

90  {
91  if(theAngularDistribution[it]!=0)
92  {
93  theAngularDistribution[it]->SetTarget(aTarget);
94  theAngularDistribution[it]->SetProjectileRP(anIncidentPart);
95  }
96  if(theEnergyAngData[it]!=0)
97  {
98  theEnergyAngData[it]->SetTarget(aTarget);
99  theEnergyAngData[it]->SetProjectileRP(anIncidentPart);
100  }
101  }
void SetProjectileRP(G4ReactionProduct &aIncidentPart)
void SetProjectileRP(const G4ReactionProduct &anIncidentParticleRP)
void SetTarget(const G4ReactionProduct &aTarget)
G4ParticleHPEnAngCorrelation * theEnergyAngData[51]
void SetTarget(G4ReactionProduct &aTarget)
G4ParticleHPAngular * theAngularDistribution[51]

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleHPInelasticCompFS::InitGammas ( G4double  AR,
G4double  ZR 
)

Definition at line 59 of file G4ParticleHPInelasticCompFS.cc.

60 {
61  // char the[100] = {""};
62  // std::ostrstream ost(the, 100, std::ios::out);
63  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
64  // G4String * aName = new G4String(the);
65  // std::ifstream from(*aName, std::ios::in);
66 
67  std::ostringstream ost;
68  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
69  G4String aName = ost.str();
70  std::ifstream from(aName, std::ios::in);
71 
72  if(!from) return; // no data found for this isotope
73  // std::ifstream theGammaData(*aName, std::ios::in);
74  std::ifstream theGammaData(aName, std::ios::in);
75 
76  theGammas.Init(theGammaData);
77  // delete aName;
78 
79 }
void Init(std::istream &aDataFile)

Here is the call graph for this function:

Here is the caller graph for this function:

G4int G4ParticleHPInelasticCompFS::SelectExitChannel ( G4double  eKinetic)

Definition at line 195 of file G4ParticleHPInelasticCompFS.cc.

196 {
197 
198 // it = 0 has without Photon
199  G4double running[50];
200  running[0] = 0;
201  unsigned int i;
202  for(i=0; i<50; i++)
203  {
204  if(i!=0) running[i]=running[i-1];
205  if(theXsection[i] != 0)
206  {
207  running[i] += std::max(0., theXsection[i]->GetXsec(eKinetic));
208  }
209  }
210  G4double random = G4UniformRand();
211  G4double sum = running[49];
212  G4int it = 50;
213  if(0!=sum)
214  {
215  G4int i0;
216  for(i0=0; i0<50; i0++)
217  {
218  it = i0;
219  // G4cout << " SelectExitChannel " << it << " " << random << " " << running[i0]/sum << " " << running[i0] << G4endl; //GDEB
220  if(random < running[i0]/sum) break;
221  }
222  }
223 //debug: it = 1;
224 // G4cout << " SelectExitChannel " << it << " " << sum << G4endl; //GDEB
225  return it;
226 }
virtual G4ParticleHPVector * GetXsec()
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

G4String G4ParticleHPInelasticCompFS::gammaPath
protected

Definition at line 113 of file G4ParticleHPInelasticCompFS.hh.

std::vector<G4int > G4ParticleHPInelasticCompFS::LR
protected

Definition at line 120 of file G4ParticleHPInelasticCompFS.hh.

std::vector< G4double > G4ParticleHPInelasticCompFS::QI
protected

Definition at line 119 of file G4ParticleHPInelasticCompFS.hh.

G4ParticleHPAngular* G4ParticleHPInelasticCompFS::theAngularDistribution[51]
protected

Definition at line 107 of file G4ParticleHPInelasticCompFS.hh.

G4ParticleHPEnAngCorrelation* G4ParticleHPInelasticCompFS::theEnergyAngData[51]
protected

Definition at line 108 of file G4ParticleHPInelasticCompFS.hh.

G4ParticleHPEnergyDistribution* G4ParticleHPInelasticCompFS::theEnergyDistribution[51]
protected

Definition at line 106 of file G4ParticleHPInelasticCompFS.hh.

G4ParticleHPPhotonDist* G4ParticleHPInelasticCompFS::theFinalStatePhotons[51]
protected

Definition at line 110 of file G4ParticleHPInelasticCompFS.hh.

G4ParticleHPDeExGammas G4ParticleHPInelasticCompFS::theGammas
protected

Definition at line 112 of file G4ParticleHPInelasticCompFS.hh.

G4ParticleHPVector* G4ParticleHPInelasticCompFS::theXsection[51]
protected

Definition at line 105 of file G4ParticleHPInelasticCompFS.hh.


The documentation for this class was generated from the following files: