Geant4  10.02.p03
G4HadronicProcess Class Reference

#include <G4HadronicProcess.hh>

Inheritance diagram for G4HadronicProcess:
Collaboration diagram for G4HadronicProcess:

Public Member Functions

 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 
virtual ~G4HadronicProcess ()
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
const G4NucleusGetTargetNucleus () const
 
const G4IsotopeGetTargetIsotope ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void BiasCrossSectionByFactor (G4double aScale)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
void MultiplyCrossSectionBy (G4double factor)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4HadronicInteractionChooseHadronicInteraction (const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, G4Material *aMaterial, G4Element *anElement)
 
G4NucleusGetTargetNucleusPointer ()
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4HadProjectile thePro
 
G4ParticleChange * theTotalResult
 
G4int epReportLevel
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Private Member Functions

G4double XBiasSurvivalProbability ()
 
G4double XBiasSecondaryWeight ()
 
G4HadronicProcessoperator= (const G4HadronicProcess &right)
 
 G4HadronicProcess (const G4HadronicProcess &)
 
void GetEnergyMomentumCheckEnvvars ()
 
G4MaterialInitialiseMaterial (G4int Z)
 

Private Attributes

G4EnergyRangeManager theEnergyRangeManager
 
G4HadronicInteractiontheInteraction
 
G4CrossSectionDataStoretheCrossSectionDataStore
 
G4HadronicProcessStoretheProcessStore
 
G4Nucleus targetNucleus
 
bool G4HadronicProcess_debug_flag
 
std::pair< G4double, G4doubleepCheckLevels
 
G4bool levelsSetByProcess
 
std::vector< G4VLeadingParticleBiasing * > theBias
 
G4double theInitialNumberOfInteractionLength
 
G4double aScaleFactor
 
G4bool xBiasOn
 
G4double theLastCrossSection
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 69 of file G4HadronicProcess.hh.

Constructor & Destructor Documentation

◆ G4HadronicProcess() [1/3]

G4HadronicProcess::G4HadronicProcess ( const G4String processName = "Hadronic",
G4ProcessType  procType = fHadronic 
)

Definition at line 86 of file G4HadronicProcess.cc.

88  : G4VDiscreteProcess(processName, procType)
89 {
90  SetProcessSubType(fHadronInelastic); // Default unless subclass changes
91 
92  theTotalResult = new G4ParticleChange();
93  theTotalResult->SetSecondaryWeightByProcess(true);
94  theInteraction = 0;
99  aScaleFactor = 1;
100  xBiasOn = false;
101  theLastCrossSection = 0.0;
104 }
static G4HadronicProcessStore * Instance()
G4HadronicInteraction * theInteraction
void GetEnergyMomentumCheckEnvvars()
G4CrossSectionDataStore * theCrossSectionDataStore
G4ParticleChange * theTotalResult
void Register(G4HadronicProcess *)
G4double theInitialNumberOfInteractionLength
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4HadronicProcessStore * theProcessStore
Here is the call graph for this function:
Here is the caller graph for this function:

◆ G4HadronicProcess() [2/3]

G4HadronicProcess::G4HadronicProcess ( const G4String processName,
G4HadronicProcessType  subType 
)

Definition at line 108 of file G4HadronicProcess.cc.

110  : G4VDiscreteProcess(processName, fHadronic)
111 {
112  SetProcessSubType(aHadSubType);
113 
114  theTotalResult = new G4ParticleChange();
115  theTotalResult->SetSecondaryWeightByProcess(true);
116  theInteraction = 0;
119  theProcessStore->Register(this);
121  aScaleFactor = 1;
122  xBiasOn = false;
123  theLastCrossSection = 0.0;
126 }
static G4HadronicProcessStore * Instance()
G4HadronicInteraction * theInteraction
void GetEnergyMomentumCheckEnvvars()
G4CrossSectionDataStore * theCrossSectionDataStore
G4ParticleChange * theTotalResult
void Register(G4HadronicProcess *)
G4double theInitialNumberOfInteractionLength
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4HadronicProcessStore * theProcessStore
Here is the call graph for this function:

◆ ~G4HadronicProcess()

G4HadronicProcess::~G4HadronicProcess ( )
virtual

Definition at line 129 of file G4HadronicProcess.cc.

130 {
132  delete theTotalResult;
134 }
G4CrossSectionDataStore * theCrossSectionDataStore
G4ParticleChange * theTotalResult
G4HadronicProcessStore * theProcessStore
void DeRegister(G4HadronicProcess *)
Here is the call graph for this function:

◆ G4HadronicProcess() [3/3]

G4HadronicProcess::G4HadronicProcess ( const G4HadronicProcess )
private

Member Function Documentation

◆ AddDataSet()

void G4HadronicProcess::AddDataSet ( G4VCrossSectionDataSet aDataSet)
inline

Definition at line 111 of file G4HadronicProcess.hh.

112  { theCrossSectionDataStore->AddDataSet(aDataSet);}
void AddDataSet(G4VCrossSectionDataSet *)
G4CrossSectionDataStore * theCrossSectionDataStore
Here is the call graph for this function:

◆ BiasCrossSectionByFactor()

void G4HadronicProcess::BiasCrossSectionByFactor ( G4double  aScale)

Definition at line 513 of file G4HadronicProcess.cc.

514 {
515  xBiasOn = true;
516  aScaleFactor = aScale;
517  G4String it = GetProcessName();
518  if ((it != "photonNuclear") &&
519  (it != "electronNuclear") &&
520  (it != "positronNuclear") ) {
522  G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had009",
523  FatalException, ed,
524  "Cross-section biasing available only for gamma and electro nuclear reactions.");
525  }
526 
527  if (aScale < 100) {
529  G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had010", JustWarning,ed,
530  "Cross-section bias readjusted to be above safe limit. New value is 100");
531  aScaleFactor = 100.;
532  }
533 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4HadronicProcess::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Reimplemented in G4HadronStoppingProcess, and G4ChargeExchangeProcess.

Definition at line 202 of file G4HadronicProcess.cc.

203 {
204  try
205  {
208  }
209  catch(G4HadronicException aR)
210  {
212  aR.Report(ed);
213  ed << " hadronic initialisation fails" << G4endl;
214  G4Exception("G4HadronicProcess::BuildPhysicsTable", "had000",
215  FatalException,ed);
216  }
218 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4HadronicProcessStore * Instance()
G4EnergyRangeManager theEnergyRangeManager
G4CrossSectionDataStore * theCrossSectionDataStore
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void BuildPhysicsTable(const G4ParticleDefinition &)
void BuildPhysicsTable(const G4ParticleDefinition &)
#define G4endl
Definition: G4ios.hh:61
void Report(std::ostream &aS)
void PrintInfo(const G4ParticleDefinition *)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CheckEnergyMomentumConservation()

void G4HadronicProcess::CheckEnergyMomentumConservation ( const G4Track &  aTrack,
const G4Nucleus aNucleus 
)
protected

Definition at line 618 of file G4HadronicProcess.cc.

620 {
621  G4int target_A=aNucleus.GetA_asInt();
622  G4int target_Z=aNucleus.GetZ_asInt();
623  G4double targetMass = G4NucleiProperties::GetNuclearMass(target_A,target_Z);
624  G4LorentzVector target4mom(0, 0, 0, targetMass);
625 
626  G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
627  G4int track_A = aTrack.GetDefinition()->GetBaryonNumber();
628  G4int track_Z = G4lrint(aTrack.GetDefinition()->GetPDGCharge());
629 
630  G4int initial_A = target_A + track_A;
631  G4int initial_Z = target_Z + track_Z;
632 
633  G4LorentzVector initial4mom = projectile4mom + target4mom;
634 
635  // Compute final-state momentum for scattering and "do nothing" results
636  G4LorentzVector final4mom;
637  G4int final_A(0), final_Z(0);
638 
639  G4int nSec = theTotalResult->GetNumberOfSecondaries();
640  if (theTotalResult->GetTrackStatus() != fStopAndKill) { // If it is Alive
641  // Either interaction didn't complete, returned "do nothing" state
642  // or the primary survived the interaction (e.g. electro-nucleus )
643  G4Track temp(aTrack);
644 
645  // Use the final energy / momentum
646  temp.SetMomentumDirection(*theTotalResult->GetMomentumDirection());
647  temp.SetKineticEnergy(theTotalResult->GetEnergy());
648 
649  if( nSec == 0 ){
650  // Interaction didn't complete, returned "do nothing" state
651  // - or suppressed recoil (e.g. Neutron elastic )
652  final4mom = temp.GetDynamicParticle()->Get4Momentum() + target4mom;
653  final_A = initial_A;
654  final_Z = initial_Z;
655  }else{
656  // The primary remains in final state (e.g. electro-nucleus )
657  final4mom = temp.GetDynamicParticle()->Get4Momentum();
658  final_A = track_A;
659  final_Z = track_Z;
660  // Expect that the target nucleus will have interacted,
661  // and its products, including recoil, will be included in secondaries.
662  }
663  }
664  if( nSec > 0 ) {
665  G4Track* sec;
666 
667  for (G4int i = 0; i < nSec; i++) {
668  sec = theTotalResult->GetSecondary(i);
669  final4mom += sec->GetDynamicParticle()->Get4Momentum();
670  final_A += sec->GetDefinition()->GetBaryonNumber();
671  final_Z += G4lrint(sec->GetDefinition()->GetPDGCharge());
672  }
673  }
674 
675  // Get level-checking information (used to cut-off relative checks)
676  G4String processName = GetProcessName();
678  G4String modelName("none");
679  if (theModel) modelName = theModel->GetModelName();
680  std::pair<G4double, G4double> checkLevels = epCheckLevels;
681  if (!levelsSetByProcess) {
682  if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
683  checkLevels.first= std::min(checkLevels.first, epCheckLevels.first);
684  checkLevels.second=std::min(checkLevels.second, epCheckLevels.second);
685  }
686 
687  // Compute absolute total-energy difference, and relative kinetic-energy
688  G4bool checkRelative = (aTrack.GetKineticEnergy() > checkLevels.second);
689 
690  G4LorentzVector diff = initial4mom - final4mom;
691  G4double absolute = diff.e();
692  G4double relative = checkRelative ? absolute/aTrack.GetKineticEnergy() : 0.;
693 
694  G4double absolute_mom = diff.vect().mag();
695  G4double relative_mom = checkRelative ? absolute_mom/aTrack.GetMomentum().mag() : 0.;
696 
697  // Evaluate relative and absolute conservation
698  G4bool relPass = true;
699  G4String relResult = "pass";
700  if ( std::abs(relative) > checkLevels.first
701  || std::abs(relative_mom) > checkLevels.first) {
702  relPass = false;
703  relResult = checkRelative ? "fail" : "N/A";
704  }
705 
706  G4bool absPass = true;
707  G4String absResult = "pass";
708  if ( std::abs(absolute) > checkLevels.second
709  || std::abs(absolute_mom) > checkLevels.second ) {
710  absPass = false ;
711  absResult = "fail";
712  }
713 
714  G4bool chargePass = true;
715  G4String chargeResult = "pass";
716  if ( (initial_A-final_A)!=0
717  || (initial_Z-final_Z)!=0 ) {
718  chargePass = checkLevels.second < DBL_MAX ? false : true;
719  chargeResult = "fail";
720  }
721 
722  G4bool conservationPass = (relPass || absPass) && chargePass;
723 
724  std::stringstream Myout;
725  G4bool Myout_notempty(false);
726  // Options for level of reporting detail:
727  // 0. off
728  // 1. report only when E/p not conserved
729  // 2. report regardless of E/p conservation
730  // 3. report only when E/p not conserved, with model names, process names, and limits
731  // 4. report regardless of E/p conservation, with model names, process names, and limits
732  // negative -1.., as above, but send output to stderr
733 
734  if( std::abs(epReportLevel) == 4
735  || ( std::abs(epReportLevel) == 3 && ! conservationPass ) ){
736  Myout << " Process: " << processName << " , Model: " << modelName << G4endl;
737  Myout << " Primary: " << aTrack.GetParticleDefinition()->GetParticleName()
738  << " (" << aTrack.GetParticleDefinition()->GetPDGEncoding() << "),"
739  << " E= " << aTrack.GetDynamicParticle()->Get4Momentum().e()
740  << ", target nucleus (" << aNucleus.GetZ_asInt() << ","
741  << aNucleus.GetA_asInt() << ")" << G4endl;
742  Myout_notempty=true;
743  }
744  if ( std::abs(epReportLevel) == 4
745  || std::abs(epReportLevel) == 2
746  || ! conservationPass ){
747 
748  Myout << " "<< relResult <<" relative, limit " << checkLevels.first << ", values E/T(0) = "
749  << relative << " p/p(0)= " << relative_mom << G4endl;
750  Myout << " "<< absResult << " absolute, limit (MeV) " << checkLevels.second/MeV << ", values E / p (MeV) = "
751  << absolute/MeV << " / " << absolute_mom/MeV << " 3mom: " << (diff.vect())*1./MeV << G4endl;
752  Myout << " "<< chargeResult << " charge/baryon number balance " << (initial_Z-final_Z) << " / " << (initial_A-final_A) << " "<< G4endl;
753  Myout_notempty=true;
754 
755  }
756  Myout.flush();
757  if ( Myout_notempty ) {
758  if (epReportLevel > 0) G4cout << Myout.str()<< G4endl;
759  else if (epReportLevel < 0) G4cerr << Myout.str()<< G4endl;
760  }
761 }
static const double MeV
Definition: G4SIunits.hh:211
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int first(char) const
G4HadronicInteraction * GetHadronicInteraction() const
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
int G4int
Definition: G4Types.hh:78
std::pair< G4double, G4double > epCheckLevels
Hep3Vector vect() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4ParticleChange * theTotalResult
bool G4bool
Definition: G4Types.hh:79
double mag() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4String & GetModelName() const
#define DBL_MAX
Definition: templates.hh:83
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CheckResult()

G4HadFinalState * G4HadronicProcess::CheckResult ( const G4HadProjectile thePro,
const G4Nucleus targetNucleus,
G4HadFinalState result 
)
protected

Definition at line 535 of file G4HadronicProcess.cc.

538 {
539  // check for catastrophic energy non-conservation
540  // to re-sample the interaction
541 
543  G4double nuclearMass(0);
544  if (theModel) {
545 
546  // Compute final-state total energy
547  G4double finalE(0.);
548  G4int nSec = result->GetNumberOfSecondaries();
549 
550  nuclearMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetA_asInt(),
551  aNucleus.GetZ_asInt());
552  if (result->GetStatusChange() != stopAndKill) {
553  // Interaction didn't complete, returned "do nothing" state
554  // and reset nucleus or the primary survived the interaction
555  // (e.g. electro-nuclear ) => keep nucleus
556  finalE=result->GetLocalEnergyDeposit() +
557  aPro.GetDefinition()->GetPDGMass() + result->GetEnergyChange();
558  if( nSec == 0 ){
559  // Since there are no secondaries, there is no recoil nucleus.
560  // To check energy balance we must neglect the initial nucleus too.
561  nuclearMass=0.0;
562  }
563  }
564  for (G4int i = 0; i < nSec; i++) {
565  G4DynamicParticle *pdyn=result->GetSecondary(i)->GetParticle();
566  finalE += pdyn->GetTotalEnergy();
567  G4double mass_pdg=pdyn->GetDefinition()->GetPDGMass();
568  G4double mass_dyn=pdyn->GetMass();
569  if ( std::abs(mass_pdg - mass_dyn) > 0.1*mass_pdg + 1.*MeV){
570  result->Clear();
571  result = 0;
573  desc << "Warning: Secondary with off-shell dynamic mass detected: " << G4endl
574  << " " << pdyn->GetDefinition()->GetParticleName()
575  << ", PDG mass: " << mass_pdg << ", dynamic mass: "<< mass_dyn << G4endl
576  << (epReportLevel<0 ? "abort the event" : "re-sample the interaction") << G4endl
577  << " Process / Model: " << GetProcessName()<< " / "
578  << theModel->GetModelName() << G4endl
579  << " Primary: " << aPro.GetDefinition()->GetParticleName()
580  << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
581  << " E= " << aPro.Get4Momentum().e()
582  << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
583  << aNucleus.GetA_asInt() << ")" << G4endl;
584  G4Exception("G4HadronicProcess:CheckResult()", "had012",
586  // must return here.....
587  return result;
588  }
589  }
590  G4double deltaE= nuclearMass + aPro.GetTotalEnergy() - finalE;
591 
592  std::pair<G4double, G4double> checkLevels =
593  theModel->GetFatalEnergyCheckLevels(); // (relative, absolute)
594  if (std::abs(deltaE) > checkLevels.second &&
595  std::abs(deltaE) > checkLevels.first*aPro.GetKineticEnergy()){
596  // do not delete result, this is a pointer to a data member;
597  result->Clear();
598  result = 0;
600  desc << "Warning: Bad energy non-conservation detected, will "
601  << (epReportLevel<0 ? "abort the event" : "re-sample the interaction") << G4endl
602  << " Process / Model: " << GetProcessName()<< " / "
603  << theModel->GetModelName() << G4endl
604  << " Primary: " << aPro.GetDefinition()->GetParticleName()
605  << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
606  << " E= " << aPro.Get4Momentum().e()
607  << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
608  << aNucleus.GetA_asInt() << ")" << G4endl
609  << " E(initial - final) = " << deltaE << " MeV." << G4endl;
610  G4Exception("G4HadronicProcess:CheckResult()", "had012",
612  }
613  }
614  return result;
615 }
G4double GetMass() const
static const double MeV
Definition: G4SIunits.hh:211
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetNumberOfSecondaries() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4HadSecondary * GetSecondary(size_t i)
G4HadronicInteraction * GetHadronicInteraction() const
G4double GetEnergyChange() const
int G4int
Definition: G4Types.hh:78
G4double GetTotalEnergy() const
G4HadFinalStateStatus GetStatusChange() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4String & GetParticleName() const
G4double GetLocalEnergyDeposit() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4DynamicParticle * GetParticle()
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4String & GetModelName() const
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ChooseHadronicInteraction()

G4HadronicInteraction* G4HadronicProcess::ChooseHadronicInteraction ( const G4HadProjectile aHadProjectile,
G4Nucleus aTargetNucleus,
G4Material aMaterial,
G4Element anElement 
)
inlineprotected

Definition at line 136 of file G4HadronicProcess.hh.

139  { return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile,
140  aTargetNucleus,
141  aMaterial,anElement);
142  }
G4EnergyRangeManager theEnergyRangeManager
G4HadronicInteraction * GetHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DumpPhysicsTable()

void G4HadronicProcess::DumpPhysicsTable ( const G4ParticleDefinition p)
inline

Definition at line 107 of file G4HadronicProcess.hh.

G4CrossSectionDataStore * theCrossSectionDataStore
void DumpPhysicsTable(const G4ParticleDefinition &)
Here is the call graph for this function:

◆ DumpState()

void G4HadronicProcess::DumpState ( const G4Track &  aTrack,
const G4String method,
G4ExceptionDescription ed 
)
protected

Definition at line 764 of file G4HadronicProcess.cc.

767 {
768  ed << "Unrecoverable error in the method " << method << " of "
769  << GetProcessName() << G4endl;
770  ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
771  << aTrack.GetParentID()
772  << " " << aTrack.GetParticleDefinition()->GetParticleName()
773  << G4endl;
774  ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
775  << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
776  ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
777 
778  if (aTrack.GetMaterial()) {
779  ed << " material " << aTrack.GetMaterial()->GetName();
780  }
781  ed << G4endl;
782 
783  if (aTrack.GetVolume()) {
784  ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
785  << ">" << G4endl;
786  }
787 }
static const double GeV
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static const double mm
Definition: SystemOfUnits.h:94
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ FillResult()

void G4HadronicProcess::FillResult ( G4HadFinalState aR,
const G4Track &  aT 
)
protected

Definition at line 412 of file G4HadronicProcess.cc.

413 {
414  theTotalResult->ProposeLocalEnergyDeposit(aR->GetLocalEnergyDeposit());
415 
416  G4double rotation = CLHEP::twopi*G4UniformRand();
417  G4ThreeVector it(0., 0., 1.);
418 
419  G4double efinal = aR->GetEnergyChange();
420  if(efinal < 0.0) { efinal = 0.0; }
421 
422  // check status of primary
423  if(aR->GetStatusChange() == stopAndKill) {
424  theTotalResult->ProposeTrackStatus(fStopAndKill);
425  theTotalResult->ProposeEnergy( 0.0 );
426 
427  // check its final energy
428  } else if(0.0 == efinal) {
429  theTotalResult->ProposeEnergy( 0.0 );
430  if(aT.GetParticleDefinition()->GetProcessManager()
431  ->GetAtRestProcessVector()->size() > 0)
432  { theTotalResult->ProposeTrackStatus(fStopButAlive); }
433  else { theTotalResult->ProposeTrackStatus(fStopAndKill); }
434 
435  // primary is not killed apply rotation and Lorentz transformation
436  } else {
437  theTotalResult->ProposeTrackStatus(fAlive);
438  G4double mass = aT.GetParticleDefinition()->GetPDGMass();
439  G4double newE = efinal + mass;
440  G4double newP = std::sqrt(efinal*(efinal + 2*mass));
441  G4ThreeVector newPV = newP*aR->GetMomentumChange();
442  G4LorentzVector newP4(newE, newPV);
443  newP4.rotate(rotation, it);
444  newP4 *= aR->GetTrafoToLab();
445  theTotalResult->ProposeMomentumDirection(newP4.vect().unit());
446  newE = newP4.e() - mass;
447  if(G4HadronicProcess_debug_flag && newE <= 0.0) {
449  DumpState(aT,"Primary has zero energy after interaction",ed);
450  G4Exception("G4HadronicProcess::FillResults", "had011", JustWarning, ed);
451  }
452  if(newE < 0.0) { newE = 0.0; }
453  theTotalResult->ProposeEnergy( newE );
454  }
455  //G4cout << "FillResult: Efinal= " << efinal << " status= "
456  // << theTotalResult->GetTrackStatus()
457  // << " fKill= " << fStopAndKill << G4endl;
458 
459  // check secondaries: apply rotation and Lorentz transformation
460  G4int nSec = aR->GetNumberOfSecondaries();
461  theTotalResult->SetNumberOfSecondaries(nSec);
462  G4double weight = aT.GetWeight();
463 
464  if (nSec > 0) {
465  G4double time0 = aT.GetGlobalTime();
466  for (G4int i = 0; i < nSec; ++i) {
468  theM.rotate(rotation, it);
469  theM *= aR->GetTrafoToLab();
470  aR->GetSecondary(i)->GetParticle()->Set4Momentum(theM);
471 
472  // time of interaction starts from zero
473  G4double time = aR->GetSecondary(i)->GetTime();
474  if (time < 0.0) { time = 0.0; }
475 
476  // take into account global time
477  time += time0;
478 
479  G4Track* track = new G4Track(aR->GetSecondary(i)->GetParticle(),
480  time, aT.GetPosition());
481  track->SetCreatorModelIndex(aR->GetSecondary(i)->GetCreatorModelType());
482  G4double newWeight = weight*aR->GetSecondary(i)->GetWeight();
483  // G4cout << "#### ParticleDebug "
484  // <<GetProcessName()<<" "
485  //<<aR->GetSecondary(i)->GetParticle()->GetDefinition()->GetParticleName()<<" "
486  // <<aScaleFactor<<" "
487  // <<XBiasSurvivalProbability()<<" "
488  // <<XBiasSecondaryWeight()<<" "
489  // <<aT.GetWeight()<<" "
490  // <<aR->GetSecondary(i)->GetWeight()<<" "
491  // <<aR->GetSecondary(i)->GetParticle()->Get4Momentum()<<" "
492  // <<G4endl;
493  track->SetWeight(newWeight);
494  track->SetTouchableHandle(aT.GetTouchableHandle());
495  theTotalResult->AddSecondary(track);
497  G4double e = track->GetKineticEnergy();
498  if (e <= 0.0) {
500  DumpState(aT,"Secondary has zero energy",ed);
501  ed << "Secondary " << track->GetDefinition()->GetParticleName()
502  << G4endl;
503  G4Exception("G4HadronicProcess::FillResults", "had011",
504  JustWarning,ed);
505  }
506  }
507  }
508  }
509 
510  aR->Clear();
511 }
G4int GetNumberOfSecondaries() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4HadSecondary * GetSecondary(size_t i)
const G4LorentzRotation & GetTrafoToLab() const
G4double GetEnergyChange() const
G4double GetWeight() const
double weight
Definition: plottest35.C:25
int G4int
Definition: G4Types.hh:78
G4HadFinalStateStatus GetStatusChange() const
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetLocalEnergyDeposit() const
G4ParticleChange * theTotalResult
G4double GetTime() const
const G4ThreeVector & GetMomentumChange() const
HepLorentzVector & rotate(double, const Hep3Vector &)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void Set4Momentum(const G4LorentzVector &momentum)
G4DynamicParticle * GetParticle()
G4LorentzVector Get4Momentum() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const double twopi
Definition: SystemOfUnits.h:54
G4int GetCreatorModelType() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCrossSectionDataStore()

G4CrossSectionDataStore* G4HadronicProcess::GetCrossSectionDataStore ( )
inline

Definition at line 166 of file G4HadronicProcess.hh.

167  {return theCrossSectionDataStore;}
G4CrossSectionDataStore * theCrossSectionDataStore
Here is the caller graph for this function:

◆ GetElementCrossSection()

G4double G4HadronicProcess::GetElementCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat = 0 
)

Definition at line 165 of file G4HadronicProcess.cc.

168 {
169  G4Material* aMaterial = const_cast<G4Material*>(mat);
170  if(! mat)
171  {
172  // Because NeutronHP needs a material pointer (for instance to get the
173  // temperature), we ask the Nist manager to build or find a simple material
174  // from the (integer) Z of the element.
175  // Note that repeated calls to this method are not producing multiple copies
176  // of the same material. But it needs to be protected against race conditions
177  // between different threads.
178  aMaterial = InitialiseMaterial(G4int(elm->GetZ()));
179  }
180  G4double x = theCrossSectionDataStore->GetCrossSection(part, elm, aMaterial);
181  if(x < 0.0) { x = 0.0; }
182  return x;
183 }
int G4int
Definition: G4Types.hh:78
Float_t mat
G4CrossSectionDataStore * theCrossSectionDataStore
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4Material * InitialiseMaterial(G4int Z)
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetEnergyMomentumCheckEnvvars()

void G4HadronicProcess::GetEnergyMomentumCheckEnvvars ( )
private

Definition at line 136 of file G4HadronicProcess.cc.

136  {
137  levelsSetByProcess = false;
138 
139  epReportLevel = getenv("G4Hadronic_epReportLevel") ?
140  strtol(getenv("G4Hadronic_epReportLevel"),0,10) : 0;
141 
142  epCheckLevels.first = getenv("G4Hadronic_epCheckRelativeLevel") ?
143  strtod(getenv("G4Hadronic_epCheckRelativeLevel"),0) : DBL_MAX;
144 
145  epCheckLevels.second = getenv("G4Hadronic_epCheckAbsoluteLevel") ?
146  strtod(getenv("G4Hadronic_epCheckAbsoluteLevel"),0) : DBL_MAX;
147 }
std::pair< G4double, G4double > epCheckLevels
#define DBL_MAX
Definition: templates.hh:83
Here is the caller graph for this function:

◆ GetEnergyMomentumCheckLevels()

std::pair<G4double, G4double> G4HadronicProcess::GetEnergyMomentumCheckLevels ( ) const
inline

Definition at line 162 of file G4HadronicProcess.hh.

163  { return epCheckLevels; }
std::pair< G4double, G4double > epCheckLevels
Here is the caller graph for this function:

◆ GetHadronicInteraction()

G4HadronicInteraction* G4HadronicProcess::GetHadronicInteraction ( ) const
inlineprotected

Definition at line 177 of file G4HadronicProcess.hh.

178  { return theInteraction; }
G4HadronicInteraction * theInteraction
Here is the caller graph for this function:

◆ GetHadronicInteractionList()

std::vector<G4HadronicInteraction*>& G4HadronicProcess::GetHadronicInteractionList ( )
inline

Definition at line 115 of file G4HadronicProcess.hh.

G4EnergyRangeManager theEnergyRangeManager
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
Here is the call graph for this function:

◆ GetLastCrossSection()

G4double G4HadronicProcess::GetLastCrossSection ( )
inlineprotected

Definition at line 181 of file G4HadronicProcess.hh.

182  { return theLastCrossSection; }
Here is the call graph for this function:

◆ GetMeanFreePath()

G4double G4HadronicProcess::GetMeanFreePath ( const G4Track &  aTrack,
G4double  ,
G4ForceCondition *   
)
virtual

Implements G4VDiscreteProcess.

Definition at line 221 of file G4HadronicProcess.cc.

222 {
223  //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
224  // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
225  try
226  {
228  theCrossSectionDataStore->GetCrossSection(aTrack.GetDynamicParticle(),
229  aTrack.GetMaterial());
230  }
231  catch(G4HadronicException aR)
232  {
234  aR.Report(ed);
235  DumpState(aTrack,"GetMeanFreePath",ed);
236  ed << " Cross section is not available" << G4endl;
237  G4Exception("G4HadronicProcess::GetMeanFreePath", "had002", FatalException,
238  ed);
239  }
240  G4double res = DBL_MAX;
241  if( theLastCrossSection > 0.0 ) { res = 1.0/theLastCrossSection; }
242  //G4cout << " xsection= " << res << G4endl;
243  return res;
244 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4CrossSectionDataStore * theCrossSectionDataStore
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void Report(std::ostream &aS)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMicroscopicCrossSection()

G4double G4HadronicProcess::GetMicroscopicCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat = 0 
)
inline

Definition at line 91 of file G4HadronicProcess.hh.

94  { return GetElementCrossSection(part, elm, mat); }
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
Here is the call graph for this function:

◆ GetTargetIsotope()

const G4Isotope* G4HadronicProcess::GetTargetIsotope ( )
inline

Definition at line 127 of file G4HadronicProcess.hh.

128  { return targetNucleus.GetIsotope(); }
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:119
Here is the call graph for this function:

◆ GetTargetNucleus()

const G4Nucleus* G4HadronicProcess::GetTargetNucleus ( ) const
inline

Definition at line 123 of file G4HadronicProcess.hh.

124  { return &targetNucleus; }

◆ GetTargetNucleusPointer()

G4Nucleus* G4HadronicProcess::GetTargetNucleusPointer ( )
inlineprotected

Definition at line 145 of file G4HadronicProcess.hh.

146  { return &targetNucleus; }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseMaterial()

G4Material * G4HadronicProcess::InitialiseMaterial ( G4int  Z)
private

Definition at line 187 of file G4HadronicProcess.cc.

188 {
189  G4AutoLock l(&hadronicProcessMutex);
191 }
G4Material * FindOrBuildSimpleMaterial(G4int Z, G4bool warning=false)
static G4NistManager * Instance()
Float_t Z
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MultiplyCrossSectionBy()

void G4HadronicProcess::MultiplyCrossSectionBy ( G4double  factor)
inline

Definition at line 169 of file G4HadronicProcess.hh.

170  { aScaleFactor = factor; }
static const G4double factor
Here is the call graph for this function:

◆ operator=()

G4HadronicProcess& G4HadronicProcess::operator= ( const G4HadronicProcess right)
private
Here is the caller graph for this function:

◆ PostStepDoIt()

G4VParticleChange * G4HadronicProcess::PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Reimplemented in CexmcHadronicProcess, and G4HadronElasticProcess.

Definition at line 247 of file G4HadronicProcess.cc.

248 {
249  //G4cout << "PostStepDoIt " << aTrack.GetDefinition()->GetParticleName()
250  // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
251  // if primary is not Alive then do nothing
252  theTotalResult->Clear();
253  theTotalResult->Initialize(aTrack);
254  theTotalResult->ProposeWeight(aTrack.GetWeight());
255  if(aTrack.GetTrackStatus() != fAlive) { return theTotalResult; }
256 
257  // Find cross section at end of step and check if <= 0
258  //
259  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
260  G4Material* aMaterial = aTrack.GetMaterial();
261 
262  G4Element* anElement = 0;
263  try
264  {
265  anElement = theCrossSectionDataStore->SampleZandA(aParticle,
266  aMaterial,
267  targetNucleus);
268  }
269  catch(G4HadronicException & aR)
270  {
272  aR.Report(ed);
273  DumpState(aTrack,"SampleZandA",ed);
274  ed << " PostStepDoIt failed on element selection" << G4endl;
275  G4Exception("G4HadronicProcess::PostStepDoIt", "had003", FatalException,
276  ed);
277  }
278 
279  // check only for charged particles
280  if(aParticle->GetDefinition()->GetPDGCharge() != 0.0) {
281  if (GetElementCrossSection(aParticle, anElement, aMaterial) <= 0.0) {
282  // No interaction
283  return theTotalResult;
284  }
285  }
286 
287  // Next check for illegal track status
288  //
289  if (aTrack.GetTrackStatus() != fAlive &&
290  aTrack.GetTrackStatus() != fSuspend) {
291  if (aTrack.GetTrackStatus() == fStopAndKill ||
292  aTrack.GetTrackStatus() == fKillTrackAndSecondaries ||
293  aTrack.GetTrackStatus() == fPostponeToNextEvent) {
295  ed << "G4HadronicProcess: track in unusable state - "
296  << aTrack.GetTrackStatus() << G4endl;
297  ed << "G4HadronicProcess: returning unchanged track " << G4endl;
298  DumpState(aTrack,"PostStepDoIt",ed);
299  G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
300  }
301  // No warning for fStopButAlive which is a legal status here
302  return theTotalResult;
303  }
304 
305  // Initialize the hadronic projectile from the track
306  thePro.Initialise(aTrack);
307 
308  try
309  {
311  ChooseHadronicInteraction( thePro, targetNucleus, aMaterial, anElement );
312  }
313  catch(G4HadronicException & aE)
314  {
316  aE.Report(ed);
317  ed << "Target element "<<anElement->GetName()<<" Z= "
318  << targetNucleus.GetZ_asInt() << " A= "
320  DumpState(aTrack,"ChooseHadronicInteraction",ed);
321  ed << " No HadronicInteraction found out" << G4endl;
322  G4Exception("G4HadronicProcess::PostStepDoIt", "had005", FatalException,
323  ed);
324  }
325 
326  G4HadFinalState* result = 0;
327  G4int reentryCount = 0;
328 
329  do
330  {
331  try
332  {
333  // Save random engine if requested for debugging
336  }
337  // Call the interaction
339  ++reentryCount;
340  }
341  catch(G4HadronicException aR)
342  {
344  aR.Report(ed);
345  ed << "Call for " << theInteraction->GetModelName() << G4endl;
346  ed << "Target element "<<anElement->GetName()<<" Z= "
348  << " A= " << targetNucleus.GetA_asInt() << G4endl;
349  DumpState(aTrack,"ApplyYourself",ed);
350  ed << " ApplyYourself failed" << G4endl;
351  G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
352  ed);
353  }
354 
355  // Check the result for catastrophic energy non-conservation
356  CheckResult(thePro, targetNucleus, result);
357 
358  if(reentryCount>100) {
360  ed << "Call for " << theInteraction->GetModelName() << G4endl;
361  ed << "Target element "<<anElement->GetName()<<" Z= "
363  << " A= " << targetNucleus.GetA_asInt() << G4endl;
364  DumpState(aTrack,"ApplyYourself",ed);
365  ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
366  G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
367  ed);
368  }
369  }
370  while(!result); /* Loop checking, 30-Oct-2015, G.Folger */
371 
373 
375 
376  FillResult(result, aTrack);
377 
378  if (epReportLevel != 0) {
380  }
381  //G4cout << "PostStepDoIt done " << G4endl;
382  return theTotalResult;
383 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4LorentzRotation & GetTrafoToLab()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
G4HadronicInteraction * theInteraction
G4HadProjectile thePro
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
int G4int
Definition: G4Types.hh:78
G4CrossSectionDataStore * theCrossSectionDataStore
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void SetTrafoToLab(const G4LorentzRotation &aT)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4ParticleChange * theTotalResult
static const char * G4Hadronic_Random_File
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void Initialise(const G4Track &aT)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
Definition: G4Element.hh:127
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
static void saveEngineStatus(const char filename[]="Config.conf")
Definition: Random.cc:176
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, G4Material *aMaterial, G4Element *anElement)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
const G4String & GetModelName() const
G4double GetPDGCharge() const
void Report(std::ostream &aS)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PreparePhysicsTable()

void G4HadronicProcess::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Reimplemented in G4HadronStoppingProcess, and G4HadronElasticProcess.

Definition at line 194 of file G4HadronicProcess.cc.

195 {
196  if(getenv("G4HadronicProcess_debug")) {
198  }
200 }
void RegisterParticle(G4HadronicProcess *, const G4ParticleDefinition *)
G4HadronicProcessStore * theProcessStore
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProcessDescription()

void G4HadronicProcess::ProcessDescription ( std::ostream &  outFile) const
virtual

Reimplemented in G4HadronStoppingProcess, G4HadronElasticProcess, G4MuonMinusCapture, G4HadronicAbsorptionFritiof, G4HadronFissionProcess, G4HadronCaptureProcess, G4PhotoNuclearProcess, G4PhotoCaptureProcess, G4HadronicAbsorptionBertini, G4PhotoFissionProcess, G4NeutronInelasticProcess, G4AntiAlphaInelasticProcess, G4AntiDeuteronInelasticProcess, G4AntiHe3InelasticProcess, G4ElectronNuclearProcess, G4MuonNuclearProcess, G4AntiNeutronInelasticProcess, G4IonInelasticProcess, G4PionMinusInelasticProcess, G4AlphaInelasticProcess, G4AntiOmegaMinusInelasticProcess, G4AntiProtonInelasticProcess, G4AntiSigmaMinusInelasticProcess, G4AntiSigmaPlusInelasticProcess, G4AntiTritonInelasticProcess, G4AntiXiMinusInelasticProcess, G4AntiXiZeroInelasticProcess, G4DeuteronInelasticProcess, G4KaonMinusInelasticProcess, G4KaonPlusInelasticProcess, G4KaonZeroLInelasticProcess, G4KaonZeroSInelasticProcess, G4LambdaInelasticProcess, G4OmegaMinusInelasticProcess, G4PionPlusInelasticProcess, G4ProtonInelasticProcess, G4SigmaMinusInelasticProcess, G4SigmaPlusInelasticProcess, G4TritonInelasticProcess, G4XiMinusInelasticProcess, G4XiZeroInelasticProcess, G4AntiLambdaInelasticProcess, G4He3InelasticProcess, and G4PositronNuclearProcess.

Definition at line 386 of file G4HadronicProcess.cc.

387 {
388  outFile << "The description for this process has not been written yet.\n";
389 }
Here is the caller graph for this function:

◆ RegisterMe()

void G4HadronicProcess::RegisterMe ( G4HadronicInteraction a)

Definition at line 149 of file G4HadronicProcess.cc.

150 {
151  if(!a) { return; }
152  try{ theEnergyRangeManager.RegisterMe( a ); }
153  catch(G4HadronicException & aE)
154  {
156  aE.Report(ed);
157  ed << "Unrecoverable error in " << GetProcessName()
158  << " to register " << a->GetModelName() << G4endl;
159  G4Exception("G4HadronicProcess::RegisterMe", "had001", FatalException,
160  ed);
161  }
163 }
void RegisterMe(G4HadronicInteraction *a)
void RegisterInteraction(G4HadronicProcess *, G4HadronicInteraction *)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4HadronicProcessStore * Instance()
G4EnergyRangeManager theEnergyRangeManager
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
const G4String & GetModelName() const
void Report(std::ostream &aS)
Here is the call graph for this function:

◆ SetEnergyMomentumCheckLevels()

void G4HadronicProcess::SetEnergyMomentumCheckLevels ( G4double  relativeLevel,
G4double  absoluteLevel 
)
inline

Definition at line 156 of file G4HadronicProcess.hh.

157  { epCheckLevels.first = relativeLevel;
158  epCheckLevels.second = absoluteLevel;
159  levelsSetByProcess = true;
160  }
std::pair< G4double, G4double > epCheckLevels
Here is the caller graph for this function:

◆ SetEpReportLevel()

void G4HadronicProcess::SetEpReportLevel ( G4int  level)
inline

Definition at line 153 of file G4HadronicProcess.hh.

154  { epReportLevel = level; }

◆ XBiasSecondaryWeight()

G4double G4HadronicProcess::XBiasSecondaryWeight ( )
private

Definition at line 402 of file G4HadronicProcess.cc.

403 {
404  G4double result = 0;
406  result =
407  1./aScaleFactor*std::exp(-nLTraversed/aScaleFactor*(1-1./aScaleFactor));
408  return result;
409 }
double G4double
Definition: G4Types.hh:76
G4double GetTotalNumberOfInteractionLengthTraversed() const
Definition: G4VProcess.hh:458
Here is the call graph for this function:
Here is the caller graph for this function:

◆ XBiasSurvivalProbability()

G4double G4HadronicProcess::XBiasSurvivalProbability ( )
private

Definition at line 392 of file G4HadronicProcess.cc.

393 {
394  G4double result = 0;
396  G4double biasedProbability = 1.-std::exp(-nLTraversed);
397  G4double realProbability = 1-std::exp(-nLTraversed/aScaleFactor);
398  result = (biasedProbability-realProbability)/biasedProbability;
399  return result;
400 }
double G4double
Definition: G4Types.hh:76
G4double GetTotalNumberOfInteractionLengthTraversed() const
Definition: G4VProcess.hh:458
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ aScaleFactor

G4double G4HadronicProcess::aScaleFactor
private

Definition at line 239 of file G4HadronicProcess.hh.

◆ epCheckLevels

std::pair<G4double, G4double> G4HadronicProcess::epCheckLevels
private

Definition at line 232 of file G4HadronicProcess.hh.

◆ epReportLevel

G4int G4HadronicProcess::epReportLevel
protected

Definition at line 215 of file G4HadronicProcess.hh.

◆ G4HadronicProcess_debug_flag

bool G4HadronicProcess::G4HadronicProcess_debug_flag
private

Definition at line 229 of file G4HadronicProcess.hh.

◆ levelsSetByProcess

G4bool G4HadronicProcess::levelsSetByProcess
private

Definition at line 233 of file G4HadronicProcess.hh.

◆ targetNucleus

G4Nucleus G4HadronicProcess::targetNucleus
private

Definition at line 227 of file G4HadronicProcess.hh.

◆ theBias

std::vector<G4VLeadingParticleBiasing *> G4HadronicProcess::theBias
private

Definition at line 235 of file G4HadronicProcess.hh.

◆ theCrossSectionDataStore

G4CrossSectionDataStore* G4HadronicProcess::theCrossSectionDataStore
private

Definition at line 223 of file G4HadronicProcess.hh.

◆ theEnergyRangeManager

G4EnergyRangeManager G4HadronicProcess::theEnergyRangeManager
private

Definition at line 219 of file G4HadronicProcess.hh.

◆ theInitialNumberOfInteractionLength

G4double G4HadronicProcess::theInitialNumberOfInteractionLength
private

Definition at line 237 of file G4HadronicProcess.hh.

◆ theInteraction

G4HadronicInteraction* G4HadronicProcess::theInteraction
private

Definition at line 221 of file G4HadronicProcess.hh.

◆ theLastCrossSection

G4double G4HadronicProcess::theLastCrossSection
private

Definition at line 241 of file G4HadronicProcess.hh.

◆ thePro

G4HadProjectile G4HadronicProcess::thePro
protected

Definition at line 211 of file G4HadronicProcess.hh.

◆ theProcessStore

G4HadronicProcessStore* G4HadronicProcess::theProcessStore
private

Definition at line 225 of file G4HadronicProcess.hh.

◆ theTotalResult

G4ParticleChange* G4HadronicProcess::theTotalResult
protected

Definition at line 213 of file G4HadronicProcess.hh.

◆ xBiasOn

G4bool G4HadronicProcess::xBiasOn
private

Definition at line 240 of file G4HadronicProcess.hh.


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