Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronicProcess Class Reference

#include <G4HadronicProcess.hh>

Inheritance diagram for G4HadronicProcess:

Public Member Functions

 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 ~G4HadronicProcess () override
void RegisterMe (G4HadronicInteraction *a)
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
void StartTracking (G4Track *track) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4ForceCondition *) override
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
void PreparePhysicsTable (const G4ParticleDefinition &) override
void BuildPhysicsTable (const G4ParticleDefinition &) override
void DumpPhysicsTable (const G4ParticleDefinition &p)
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
G4HadronicInteractionGetHadronicModel (const G4String &)
G4HadronicInteractionGetHadronicInteraction () const
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
const G4NucleusGetTargetNucleus () const
G4NucleusGetTargetNucleusPointer ()
const G4IsotopeGetTargetIsotope ()
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
G4HadXSType CrossSectionType () const
void SetCrossSectionType (G4HadXSType val)
void ProcessDescription (std::ostream &outFile) const override
void BiasCrossSectionByFactor (G4double aScale)
void MultiplyCrossSectionBy (G4double factor)
G4double CrossSectionFactor () const
void SetIntegral (G4bool val)
void SetEpReportLevel (G4int level)
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
G4CrossSectionDataStoreGetCrossSectionDataStore ()
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS () const
std::vector< G4double > * EnergyOfCrossSectionMax () const
G4HadronicProcessoperator= (const G4HadronicProcess &right)=delete
 G4HadronicProcess (const G4HadronicProcess &)=delete
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 G4VProcess (const G4VProcess &right)
virtual ~G4VProcess ()
G4VProcessoperator= (const G4VProcess &)=delete
G4bool operator== (const G4VProcess &right) const
G4bool 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 const G4VProcessGetCreatorProcess () const
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, const G4Material *aMaterial, const G4Element *anElement)
G4double GetLastCrossSection ()
void FillResult (G4HadFinalState *aR, const G4Track &aT)
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Protected Attributes

G4HadProjectile thePro
G4ParticleChangetheTotalResult
G4CrossSectionDataStoretheCrossSectionDataStore
G4double fWeight = 1.0
G4double aScaleFactor = 1.0
G4double theLastCrossSection = 0.0
G4double mfpKinEnergy = DBL_MAX
G4int epReportLevel = 0
G4HadXSType fXSType = fHadNoIntegral
Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
G4VParticleChangepParticleChange = nullptr
G4ParticleChange aParticleChange
G4double theNumberOfInteractionLengthLeft = -1.0
G4double currentInteractionLength = -1.0
G4double theInitialNumberOfInteractionLength = -1.0
G4String theProcessName
G4String thePhysicsTableFileName
G4ProcessType theProcessType = fNotDefined
G4int theProcessSubType = -1
G4double thePILfactor = 1.0
G4int verboseLevel = 0
G4bool enableAtRestDoIt = true
G4bool enableAlongStepDoIt = true
G4bool enablePostStepDoIt = true

Additional Inherited Members

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

Detailed Description

Definition at line 73 of file G4HadronicProcess.hh.

Constructor & Destructor Documentation

◆ G4HadronicProcess() [1/3]

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

◆ G4HadronicProcess() [2/3]

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

Definition at line 95 of file G4HadronicProcess.cc.

97 : G4VDiscreteProcess(processName, fHadronic)
98{
99 SetProcessSubType(aHadSubType);
100 InitialiseLocal();
101}
@ fHadronic

◆ ~G4HadronicProcess()

G4HadronicProcess::~G4HadronicProcess ( )
override

Definition at line 103 of file G4HadronicProcess.cc.

104{
105 theProcessStore->DeRegister(this);
106 delete theTotalResult;
108 if(isMaster) {
109 if (fXSpeaks != nullptr) {
110 for (auto const& e : *fXSpeaks ) {
111 delete e;
112 }
113 }
114 delete fXSpeaks;
115 delete theEnergyOfCrossSectionMax;
116 }
117}
G4ParticleChange * theTotalResult
G4CrossSectionDataStore * theCrossSectionDataStore

◆ G4HadronicProcess() [3/3]

G4HadronicProcess::G4HadronicProcess ( const G4HadronicProcess & )
delete

Member Function Documentation

◆ AddDataSet()

void G4HadronicProcess::AddDataSet ( G4VCrossSectionDataSet * aDataSet)

Definition at line 854 of file G4HadronicProcess.cc.

855{
856 theCrossSectionDataStore->AddDataSet(aDataSet);
857}

Referenced by G4HadProcesses::AddCaptureCrossSection(), G4HadProcesses::AddElasticCrossSection(), G4HadProcesses::AddFissionCrossSection(), G4HadProcesses::AddInelasticCrossSection(), G4HadronElasticPhysics::AddXSection(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPiKBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4NeutronPHPBuilder::Build(), G4NeutronPHPBuilder::Build(), G4NeutronPHPBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4QGSPProtonBuilder::Build(), G4TritonPHPBuilder::Build(), G4HadProcesses::BuildNeutronElastic(), G4HadProcesses::BuildNeutronInelasticAndCapture(), LBE::ConstructHad(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronElasticPhysicsHP::ConstructProcess(), G4HadronElasticPhysicsHPT::ConstructProcess(), G4HadronElasticPhysicsLEND::ConstructProcess(), G4HadronElasticPhysicsPHP::ConstructProcess(), G4HadronElasticPhysicsVI::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4HadronInelasticQBBC::ConstructProcess(), G4HadronInelasticQBBC_ABLA::ConstructProcess(), G4HadronPhysicsFTFQGSP_BERT::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4NeutronCrossSectionXS::ConstructProcess(), G4ThermalNeutrons::ConstructProcess(), G4URRNeutrons::ConstructProcess(), G4ChargeExchangeProcess::G4ChargeExchangeProcess(), G4MuonNuclearProcess::G4MuonNuclearProcess(), G4NeutronCaptureProcess::G4NeutronCaptureProcess(), G4NeutronFissionProcess::G4NeutronFissionProcess(), G4HadronPhysicsFTF_BIC::Neutron(), G4HadronPhysicsFTFP_BERT::Neutron(), G4NeutronGeneralProcess::SetCaptureProcess(), G4NeutronGeneralProcess::SetElasticProcess(), and G4NeutronGeneralProcess::SetInelasticProcess().

◆ BiasCrossSectionByFactor()

void G4HadronicProcess::BiasCrossSectionByFactor ( G4double aScale)

Definition at line 576 of file G4HadronicProcess.cc.

577{
578 if (aScale <= 0.0) {
580 ed << " Wrong biasing factor " << aScale << " for " << GetProcessName();
581 G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had010",
582 JustWarning, ed, "Cross-section bias is ignored");
583 } else {
584 aScaleFactor = aScale;
585 }
586}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
const G4String & GetProcessName() const

Referenced by MultiplyCrossSectionBy().

◆ BuildPhysicsTable()

void G4HadronicProcess::BuildPhysicsTable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4HadronStoppingProcess, and G4NeutronGeneralProcess.

Definition at line 171 of file G4HadronicProcess.cc.

172{
173 if(firstParticle != &p) { return; }
174
175 theCrossSectionDataStore->BuildPhysicsTable(p);
176 theEnergyRangeManager.BuildPhysicsTable(p);
177 G4HadronicParameters* param = G4HadronicParameters::Instance();
178
179 G4int subtype = GetProcessSubType();
180 if(useIntegralXS) {
181 if(subtype == fHadronInelastic) {
182 useIntegralXS = param->EnableIntegralInelasticXS();
183 } else if(subtype == fHadronElastic) {
184 useIntegralXS = param->EnableIntegralElasticXS();
185 }
186 }
188
189 if(nullptr == masterProcess) {
190 masterProcess = dynamic_cast<const G4HadronicProcess*>(GetMasterProcess());
191 }
192
193 // check particle for integral method
194 if(isMaster || nullptr == masterProcess) {
195 G4double charge = p.GetPDGCharge()/eplus;
196
197 // select cross section shape
198 if(charge != 0.0 && useIntegralXS) {
199 G4double tmax = param->GetMaxEnergy();
200 currentParticle = firstParticle;
201 // initialisation in the master thread
202 G4int pdg = p.GetPDGEncoding();
203 if (std::abs(pdg) == 211) {
205 } else if (pdg == 321) {
207 } else if (pdg == -321) {
209 } else if (pdg == 2212) {
211 } else if (pdg == -2212 || pdg == -1000010020 || pdg == -1000010030 ||
212 pdg == -1000020030 || pdg == -1000020040) {
214 } else if (charge > 0.0 || pdg == 11 || pdg == 13) {
216 }
217
218 delete theEnergyOfCrossSectionMax;
219 theEnergyOfCrossSectionMax = nullptr;
220 if(fXSType == fHadTwoPeaks) {
221 if (fXSpeaks != nullptr) {
222 for (auto const& e : *fXSpeaks ) {
223 delete e;
224 }
225 }
226 delete fXSpeaks;
227 fXSpeaks =
228 G4HadXSHelper::FillPeaksStructure(this, &p, minKinEnergy, tmax);
229 if(nullptr == fXSpeaks) {
231 }
232 }
233 if(fXSType == fHadOnePeak) {
234 theEnergyOfCrossSectionMax =
235 G4HadXSHelper::FindCrossSectionMax(this, &p, minKinEnergy, tmax);
236 if(nullptr == theEnergyOfCrossSectionMax) {
238 }
239 }
240 }
241 } else {
242 // initialisation in worker threads
243 fXSType = masterProcess->CrossSectionType();
244 fXSpeaks = masterProcess->TwoPeaksXS();
245 theEnergyOfCrossSectionMax = masterProcess->EnergyOfCrossSectionMax();
246 }
247 if(isMaster && 1 < param->GetVerboseLevel()) {
248 G4cout << "G4HadronicProcess::BuildPhysicsTable: for "
249 << GetProcessName() << " and " << p.GetParticleName()
250 << " typeXS=" << fXSType << G4endl;
251 }
253}
@ fHadTwoPeaks
@ fHadIncreasing
@ fHadDecreasing
@ fHadNoIntegral
@ fHadOnePeak
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static std::vector< G4TwoPeaksHadXS * > * FillPeaksStructure(G4HadronicProcess *, const G4ParticleDefinition *, const G4double tmin, const G4double tmax)
static std::vector< G4double > * FindCrossSectionMax(G4HadronicProcess *, const G4ParticleDefinition *, const G4double tmin, const G4double tmax)
G4bool EnableIntegralInelasticXS() const
static G4HadronicParameters * Instance()
G4bool EnableIntegralElasticXS() const
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
const G4String & GetParticleName() const
const G4VProcess * GetMasterProcess() const
G4int GetProcessSubType() const

Referenced by G4ChargeExchangeProcess::BuildPhysicsTable().

◆ CheckEnergyMomentumConservation()

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

Definition at line 679 of file G4HadronicProcess.cc.

681{
682 G4int target_A=aNucleus.GetA_asInt();
683 G4int target_Z=aNucleus.GetZ_asInt();
684 G4double targetMass = G4NucleiProperties::GetNuclearMass(target_A,target_Z);
685 G4LorentzVector target4mom(0, 0, 0, targetMass
686 + nICelectrons*CLHEP::electron_mass_c2);
687
688 G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
689 G4int track_A = aTrack.GetDefinition()->GetBaryonNumber();
690 G4int track_Z = G4lrint(aTrack.GetDefinition()->GetPDGCharge());
691
692 G4int initial_A = target_A + track_A;
693 G4int initial_Z = target_Z + track_Z - nICelectrons;
694
695 G4LorentzVector initial4mom = projectile4mom + target4mom;
696
697 // Compute final-state momentum for scattering and "do nothing" results
698 G4LorentzVector final4mom;
699 G4int final_A(0), final_Z(0);
700
701 G4int nSec = theTotalResult->GetNumberOfSecondaries();
702 if (theTotalResult->GetTrackStatus() != fStopAndKill) { // If it is Alive
703 // Either interaction didn't complete, returned "do nothing" state
704 // or the primary survived the interaction (e.g. electro-nucleus )
705
706 // Interaction didn't complete, returned "do nothing" state
707 // - or suppressed recoil (e.g. Neutron elastic )
708 final4mom = initial4mom;
709 final_A = initial_A;
710 final_Z = initial_Z;
711 if (nSec > 0) {
712 // The primary remains in final state (e.g. electro-nucleus )
713 // Use the final energy / momentum
714 const G4ThreeVector& v = *theTotalResult->GetMomentumDirection();
715 G4double ekin = theTotalResult->GetEnergy();
716 G4double mass = aTrack.GetDefinition()->GetPDGMass();
717 G4double ptot = std::sqrt(ekin*(ekin + 2*mass));
718 final4mom.set(ptot*v.x(), ptot*v.y(), ptot*v.z(), mass + ekin);
719 final_A = track_A;
720 final_Z = track_Z;
721 // Expect that the target nucleus will have interacted,
722 // and its products, including recoil, will be included in secondaries.
723 }
724 }
725 if( nSec > 0 ) {
726 G4Track* sec;
727
728 for (G4int i = 0; i < nSec; i++) {
729 sec = theTotalResult->GetSecondary(i);
730 final4mom += sec->GetDynamicParticle()->Get4Momentum();
731 final_A += sec->GetDefinition()->GetBaryonNumber();
732 final_Z += G4lrint(sec->GetDefinition()->GetPDGCharge());
733 }
734 }
735
736 // Get level-checking information (used to cut-off relative checks)
737 G4String processName = GetProcessName();
738 G4HadronicInteraction* theModel = GetHadronicInteraction();
739 G4String modelName("none");
740 if (theModel) modelName = theModel->GetModelName();
741 std::pair<G4double, G4double> checkLevels = epCheckLevels;
742 if (!levelsSetByProcess) {
743 if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
744 checkLevels.first= std::min(checkLevels.first, epCheckLevels.first);
745 checkLevels.second=std::min(checkLevels.second, epCheckLevels.second);
746 }
747
748 // Compute absolute total-energy difference, and relative kinetic-energy
749 G4bool checkRelative = (aTrack.GetKineticEnergy() > checkLevels.second);
750
751 G4LorentzVector diff = initial4mom - final4mom;
752 G4double absolute = diff.e();
753 G4double relative = checkRelative ? absolute/aTrack.GetKineticEnergy() : 0.;
754
755 G4double absolute_mom = diff.vect().mag();
756 G4double relative_mom = checkRelative ? absolute_mom/aTrack.GetMomentum().mag() : 0.;
757
758 // Evaluate relative and absolute conservation
759 G4bool relPass = true;
760 G4String relResult = "pass";
761 if ( std::abs(relative) > checkLevels.first
762 || std::abs(relative_mom) > checkLevels.first) {
763 relPass = false;
764 relResult = checkRelative ? "fail" : "N/A";
765 }
766
767 G4bool absPass = true;
768 G4String absResult = "pass";
769 if ( std::abs(absolute) > checkLevels.second
770 || std::abs(absolute_mom) > checkLevels.second ) {
771 absPass = false ;
772 absResult = "fail";
773 }
774
775 G4bool chargePass = true;
776 G4String chargeResult = "pass";
777 if ( (initial_A-final_A)!=0
778 || (initial_Z-final_Z)!=0 ) {
779 chargePass = checkLevels.second < DBL_MAX ? false : true;
780 chargeResult = "fail";
781 }
782
783 G4bool conservationPass = (relPass || absPass) && chargePass;
784
785 std::stringstream Myout;
786 G4bool Myout_notempty(false);
787 // Options for level of reporting detail:
788 // 0. off
789 // 1. report only when E/p not conserved
790 // 2. report regardless of E/p conservation
791 // 3. report only when E/p not conserved, with model names, process names, and limits
792 // 4. report regardless of E/p conservation, with model names, process names, and limits
793 // negative -1.., as above, but send output to stderr
794
795 if( std::abs(epReportLevel) == 4
796 || ( std::abs(epReportLevel) == 3 && ! conservationPass ) ){
797 Myout << " Process: " << processName << " , Model: " << modelName << G4endl;
798 Myout << " Primary: " << aTrack.GetParticleDefinition()->GetParticleName()
799 << " (" << aTrack.GetParticleDefinition()->GetPDGEncoding() << "),"
800 << " E= " << aTrack.GetDynamicParticle()->Get4Momentum().e()
801 << ", target nucleus (" << aNucleus.GetZ_asInt() << ","
802 << aNucleus.GetA_asInt() << ")" << G4endl;
803 Myout_notempty=true;
804 }
805 if ( std::abs(epReportLevel) == 4
806 || std::abs(epReportLevel) == 2
807 || ! conservationPass ){
808
809 Myout << " "<< relResult <<" relative, limit " << checkLevels.first << ", values E/T(0) = "
810 << relative << " p/p(0)= " << relative_mom << G4endl;
811 Myout << " "<< absResult << " absolute, limit (MeV) " << checkLevels.second/MeV << ", values E / p (MeV) = "
812 << absolute/MeV << " / " << absolute_mom/MeV << " 3mom: " << (diff.vect())*1./MeV << G4endl;
813 Myout << " "<< chargeResult << " charge/baryon number balance " << (initial_Z-final_Z) << " / " << (initial_A-final_A) << " "<< G4endl;
814 Myout_notempty=true;
815
816 }
817 Myout.flush();
818 if ( Myout_notempty ) {
819 if (epReportLevel > 0) G4cout << Myout.str()<< G4endl;
820 else if (epReportLevel < 0) G4cerr << Myout.str()<< G4endl;
821 }
822}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
bool G4bool
Definition G4Types.hh:86
G4GLOB_DLL std::ostream G4cerr
double z() const
double x() const
double y() const
double mag() const
Hep3Vector vect() const
void set(double x, double y, double z, double t)
G4LorentzVector Get4Momentum() const
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
const G4String & GetModelName() const
G4HadronicInteraction * GetHadronicInteraction() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition G4Nucleus.hh:78
G4int GetZ_asInt() const
Definition G4Nucleus.hh:84
const G4ParticleDefinition * GetParticleDefinition() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62

Referenced by G4HadronStoppingProcess::AtRestDoIt(), G4HadronElasticProcess::PostStepDoIt(), and PostStepDoIt().

◆ CheckResult()

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

Definition at line 588 of file G4HadronicProcess.cc.

591{
592 // check for catastrophic energy non-conservation
593 // to re-sample the interaction
594 G4HadronicInteraction * theModel = GetHadronicInteraction();
595 G4double nuclearMass(0);
596 if (nullptr != theModel) {
597
598 // Compute final-state total energy
599 G4double finalE(0.);
600 G4int nSec = (G4int)result->GetNumberOfSecondaries();
601
602 nuclearMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetA_asInt(),
603 aNucleus.GetZ_asInt());
604 if (result->GetStatusChange() != stopAndKill) {
605 // Interaction didn't complete, returned "do nothing" state
606 // and reset nucleus or the primary survived the interaction
607 // (e.g. electro-nuclear ) => keep nucleus
608 finalE=result->GetLocalEnergyDeposit() +
609 aPro.GetDefinition()->GetPDGMass() + result->GetEnergyChange();
610 if( nSec == 0 ){
611 // Since there are no secondaries, there is no recoil nucleus.
612 // To check energy balance we must neglect the initial nucleus too.
613 nuclearMass=0.0;
614 }
615 }
616 for (G4int i = 0; i < nSec; ++i) {
617 G4DynamicParticle *pdyn=result->GetSecondary(i)->GetParticle();
618 finalE += pdyn->GetTotalEnergy();
619 G4double mass_pdg=pdyn->GetDefinition()->GetPDGMass();
620 G4double mass_dyn=pdyn->GetMass();
621 if ( std::abs(mass_pdg - mass_dyn) > 0.1*mass_pdg + 1.*MeV ) {
622 // If it is shortlived, then a difference less than 3 times the width is acceptable
623 if ( pdyn->GetDefinition()->IsShortLived() &&
624 std::abs(mass_pdg - mass_dyn) < 3.0*pdyn->GetDefinition()->GetPDGWidth() ) {
625 continue;
626 }
627 result->Clear();
628 result = nullptr;
630 desc << "Warning: Secondary with off-shell dynamic mass detected: "
631 << G4endl
632 << " " << pdyn->GetDefinition()->GetParticleName()
633 << ", PDG mass: " << mass_pdg << ", dynamic mass: "
634 << mass_dyn << G4endl
635 << (epReportLevel<0 ? "abort the event"
636 : "re-sample the interaction") << G4endl
637 << " Process / Model: " << GetProcessName()<< " / "
638 << theModel->GetModelName() << G4endl
639 << " Primary: " << aPro.GetDefinition()->GetParticleName()
640 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
641 << " E= " << aPro.Get4Momentum().e()
642 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
643 << aNucleus.GetA_asInt() << ")" << G4endl;
644 G4Exception("G4HadronicProcess:CheckResult()", "had012",
646 // must return here.....
647 return result;
648 }
649 }
650 G4double deltaE= nuclearMass + aPro.GetTotalEnergy() - finalE;
651
652 std::pair<G4double, G4double> checkLevels =
653 theModel->GetFatalEnergyCheckLevels(); // (relative, absolute)
654 if (std::abs(deltaE) > checkLevels.second &&
655 std::abs(deltaE) > checkLevels.first*aPro.GetKineticEnergy()){
656 // do not delete result, this is a pointer to a data member;
657 result->Clear();
658 result = nullptr;
660 desc << "Warning: Bad energy non-conservation detected, will "
661 << (epReportLevel<0 ? "abort the event"
662 : "re-sample the interaction") << G4endl
663 << " Process / Model: " << GetProcessName()<< " / "
664 << theModel->GetModelName() << G4endl
665 << " Primary: " << aPro.GetDefinition()->GetParticleName()
666 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
667 << " E= " << aPro.Get4Momentum().e()
668 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
669 << aNucleus.GetA_asInt() << ")" << G4endl
670 << " E(initial - final) = " << deltaE << " MeV." << G4endl;
671 G4Exception("G4HadronicProcess:CheckResult()", "had012",
673 }
674 }
675 return result;
676}
@ EventMustBeAborted
@ stopAndKill
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4double GetEnergyChange() const
G4HadFinalStateStatus GetStatusChange() const
G4double GetLocalEnergyDeposit() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
G4DynamicParticle * GetParticle()
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const

Referenced by G4HadronStoppingProcess::AtRestDoIt(), G4HadronElasticProcess::PostStepDoIt(), and PostStepDoIt().

◆ ChooseHadronicInteraction()

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

Definition at line 349 of file G4HadronicProcess.hh.

354{
355 return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile,
356 aTargetNucleus,
357 aMaterial,anElement);
358}

Referenced by G4HadronStoppingProcess::AtRestDoIt(), G4HadronElasticProcess::PostStepDoIt(), and PostStepDoIt().

◆ ComputeCrossSection()

G4double G4HadronicProcess::ComputeCrossSection ( const G4ParticleDefinition * part,
const G4Material * mat,
const G4double kinEnergy )

Definition at line 877 of file G4HadronicProcess.cc.

880{
881 auto dp = new G4DynamicParticle(part, unitVector, kinEnergy);
882 G4double xs = theCrossSectionDataStore->ComputeCrossSection(dp, mat);
883 delete dp;
884 return xs;
885}

Referenced by G4HadXSHelper::FillPeaksStructure(), and G4HadXSHelper::FindCrossSectionMax().

◆ CrossSectionFactor()

G4double G4HadronicProcess::CrossSectionFactor ( ) const
inline

Definition at line 301 of file G4HadronicProcess.hh.

302{
303 return aScaleFactor;
304}

◆ CrossSectionType()

G4HadXSType G4HadronicProcess::CrossSectionType ( ) const
inline

Definition at line 290 of file G4HadronicProcess.hh.

291{
292 return fXSType;
293}

◆ DumpPhysicsTable()

void G4HadronicProcess::DumpPhysicsTable ( const G4ParticleDefinition & p)

Definition at line 849 of file G4HadronicProcess.cc.

850{
851 theCrossSectionDataStore->DumpPhysicsTable(p);
852}

◆ DumpState()

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

Definition at line 824 of file G4HadronicProcess.cc.

827{
828 ed << "Unrecoverable error in the method " << method << " of "
829 << GetProcessName() << G4endl;
830 ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
831 << aTrack.GetParentID()
832 << " " << aTrack.GetParticleDefinition()->GetParticleName()
833 << G4endl;
834 ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
835 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
836 ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
837
838 if (aTrack.GetMaterial()) {
839 ed << " material " << aTrack.GetMaterial()->GetName();
840 }
841 ed << G4endl;
842
843 if (aTrack.GetVolume()) {
844 ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
845 << ">" << G4endl;
846 }
847}
const G4String & GetName() const
G4int GetTrackID() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4Material * GetMaterial() const
const G4ThreeVector & GetMomentumDirection() const
G4int GetParentID() const
const G4String & GetName() const

Referenced by G4HadronStoppingProcess::AtRestDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), and G4TauNeutrinoNucleusProcess::PostStepDoIt().

◆ EnergyOfCrossSectionMax()

std::vector< G4double > * G4HadronicProcess::EnergyOfCrossSectionMax ( ) const
inline

Definition at line 344 of file G4HadronicProcess.hh.

345{
346 return theEnergyOfCrossSectionMax;
347}

◆ FillResult()

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

Definition at line 482 of file G4HadronicProcess.cc.

483{
484 theTotalResult->ProposeLocalEnergyDeposit(aR->GetLocalEnergyDeposit());
485 const G4ThreeVector& dir = aT.GetMomentumDirection();
486
487 G4double efinal = std::max(aR->GetEnergyChange(), 0.0);
488
489 // check status of primary
490 if(aR->GetStatusChange() == stopAndKill) {
491 theTotalResult->ProposeTrackStatus(fStopAndKill);
492 theTotalResult->ProposeEnergy( 0.0 );
493
494 // check its final energy
495 } else if(0.0 == efinal) {
496 theTotalResult->ProposeEnergy( 0.0 );
498 ->GetAtRestProcessVector()->size() > 0)
499 { theTotalResult->ProposeTrackStatus(fStopButAlive); }
500 else { theTotalResult->ProposeTrackStatus(fStopAndKill); }
501
502 // primary is not killed apply rotation and Lorentz transformation
503 } else {
504 theTotalResult->ProposeTrackStatus(fAlive);
505 G4ThreeVector newDir = aR->GetMomentumChange();
506 newDir.rotateUz(dir);
507 theTotalResult->ProposeMomentumDirection(newDir);
508 theTotalResult->ProposeEnergy(efinal);
509 }
510 //G4cout << "FillResult: Efinal= " << efinal << " status= "
511 // << theTotalResult->GetTrackStatus()
512 // << " fKill= " << fStopAndKill << G4endl;
513
514 // check secondaries
515 nICelectrons = 0;
516 G4int nSec = (G4int)aR->GetNumberOfSecondaries();
517 theTotalResult->SetNumberOfSecondaries(nSec);
518 G4double time0 = aT.GetGlobalTime();
519
520 for (G4int i = 0; i < nSec; ++i) {
521 G4DynamicParticle* dynParticle = aR->GetSecondary(i)->GetParticle();
522
523 // apply rotation
524 G4ThreeVector newDir = dynParticle->GetMomentumDirection();
525 newDir.rotateUz(dir);
526 dynParticle->SetMomentumDirection(newDir);
527
528 // check if secondary is on the mass shell
529 const G4ParticleDefinition* part = dynParticle->GetDefinition();
530 G4double mass = part->GetPDGMass();
531 G4double dmass= dynParticle->GetMass();
532 const G4double delta_mass_lim = 1.0*CLHEP::keV;
533 const G4double delta_ekin = 0.001*CLHEP::eV;
534 if(std::abs(dmass - mass) > delta_mass_lim) {
535 G4double e =
536 std::max(dynParticle->GetKineticEnergy() + dmass - mass, delta_ekin);
537 if(verboseLevel > 1) {
539 ed << "TrackID= "<< aT.GetTrackID()
540 << " " << aT.GetParticleDefinition()->GetParticleName()
541 << " Target Z= " << targetNucleus.GetZ_asInt() << " A= "
542 << targetNucleus.GetA_asInt()
543 << " Ekin(GeV)= " << aT.GetKineticEnergy()/CLHEP::GeV
544 << "\n Secondary is out of mass shell: " << part->GetParticleName()
545 << " EkinNew(MeV)= " << e
546 << " DeltaMass(MeV)= " << dmass - mass << G4endl;
547 G4Exception("G4HadronicProcess::FillResults", "had012", JustWarning, ed);
548 }
549 dynParticle->SetKineticEnergy(e);
550 dynParticle->SetMass(mass);
551 }
552 G4int idModel = aR->GetSecondary(i)->GetCreatorModelID();
553 if(part->GetPDGEncoding() == 11) { ++nICelectrons; }
554
555 // time of interaction starts from zero + global time
556 G4double time = std::max(aR->GetSecondary(i)->GetTime(), 0.0) + time0;
557
558 G4Track* track = new G4Track(dynParticle, time, aT.GetPosition());
559 track->SetCreatorModelID(idModel);
562 G4double newWeight = fWeight*aR->GetSecondary(i)->GetWeight();
563 track->SetWeight(newWeight);
565 theTotalResult->AddSecondary(track);
566 }
567 aR->Clear();
568 // G4cout << "FillResults done nICe= " << nICelectrons << G4endl;
569}
@ fAlive
@ fStopButAlive
Hep3Vector & rotateUz(const Hep3Vector &)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
void SetMass(G4double mass)
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
const G4ThreeVector & GetMomentumChange() const
G4int GetParentResonanceID() const
G4double GetWeight() const
const G4ParticleDefinition * GetParentResonanceDef() const
G4double GetTime() const
G4int GetCreatorModelID() const
G4ProcessManager * GetProcessManager() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
void SetWeight(G4double aValue)
void SetParentResonanceID(const G4int parentID)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
void SetParentResonanceDef(const G4ParticleDefinition *parent)
G4int verboseLevel

Referenced by G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), and G4TauNeutrinoNucleusProcess::PostStepDoIt().

◆ GetCrossSectionDataStore()

◆ GetElementCrossSection()

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

Definition at line 144 of file G4HadronicProcess.cc.

147{
148 if(nullptr == mat)
149 {
150 static const G4int nmax = 5;
151 if(nMatWarn < nmax) {
152 ++nMatWarn;
154 ed << "Cannot compute Element x-section for " << GetProcessName()
155 << " because no material defined \n"
156 << " Please, specify material pointer or define simple material"
157 << " for Z= " << elm->GetZasInt();
158 G4Exception("G4HadronicProcess::GetElementCrossSection", "had066",
159 JustWarning, ed);
160 }
161 }
162 return theCrossSectionDataStore->GetCrossSection(dp, elm, mat);
163}
G4int GetZasInt() const
Definition G4Element.hh:120

Referenced by G4HadronicProcessStore::GetCaptureCrossSectionPerAtom(), G4HadronicProcessStore::GetChargeExchangeCrossSectionPerAtom(), G4HadronicProcessStore::GetElasticCrossSectionPerAtom(), G4HadronicProcessStore::GetFissionCrossSectionPerAtom(), G4HadronicProcessStore::GetInelasticCrossSectionPerAtom(), and GetMicroscopicCrossSection().

◆ GetEnergyMomentumCheckLevels()

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

Definition at line 326 of file G4HadronicProcess.hh.

327{
328 return epCheckLevels;
329}

◆ GetHadronicInteraction()

G4HadronicInteraction * G4HadronicProcess::GetHadronicInteraction ( ) const
inline

Definition at line 273 of file G4HadronicProcess.hh.

274{
275 return theInteraction;
276}

Referenced by CheckEnergyMomentumConservation(), and CheckResult().

◆ GetHadronicInteractionList()

◆ GetHadronicModel()

G4HadronicInteraction * G4HadronicProcess::GetHadronicModel ( const G4String & modelName)

Definition at line 866 of file G4HadronicProcess.cc.

867{
868 std::vector<G4HadronicInteraction*>& list
869 = theEnergyRangeManager.GetHadronicInteractionList();
870 for (auto & mod : list) {
871 if (mod->GetModelName() == modelName) return mod;
872 }
873 return nullptr;
874}

◆ GetLastCrossSection()

G4double G4HadronicProcess::GetLastCrossSection ( )
inlineprotected

Definition at line 365 of file G4HadronicProcess.hh.

366{
367 return theLastCrossSection;
368}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Reimplemented in G4MuNeutrinoNucleusProcess, G4NeutrinoElectronProcess, G4NeutronGeneralProcess, and G4TauNeutrinoNucleusProcess.

Definition at line 299 of file G4HadronicProcess.cc.

302{
304 ->ComputeCrossSection(aTrack.GetDynamicParticle(),aTrack.GetMaterial());
305 return (xs > 0.0) ? 1.0/xs : DBL_MAX;
306}

◆ GetMicroscopicCrossSection()

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

Definition at line 264 of file G4HadronicProcess.hh.

268{
269 return GetElementCrossSection(part, elm, mat);
270}
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)

◆ GetTargetIsotope()

const G4Isotope * G4HadronicProcess::GetTargetIsotope ( )
inline

Definition at line 284 of file G4HadronicProcess.hh.

285{
286 return targetNucleus.GetIsotope();
287}

◆ GetTargetNucleus()

const G4Nucleus * G4HadronicProcess::GetTargetNucleus ( ) const
inline

Definition at line 279 of file G4HadronicProcess.hh.

280{
281 return &targetNucleus;
282}

◆ GetTargetNucleusPointer()

◆ MultiplyCrossSectionBy()

◆ operator=()

G4HadronicProcess & G4HadronicProcess::operator= ( const G4HadronicProcess & right)
delete

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Reimplemented in G4MuNeutrinoNucleusProcess, G4NeutrinoElectronProcess, G4NeutronGeneralProcess, and G4TauNeutrinoNucleusProcess.

Definition at line 309 of file G4HadronicProcess.cc.

310{
312
313 //G4cout << "PostStepDoIt " << aTrack.GetDefinition()->GetParticleName()
314 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
315 // if primary is not Alive then do nothing
316 theTotalResult->Clear();
317 theTotalResult->Initialize(aTrack);
318 fWeight = aTrack.GetWeight();
319 theTotalResult->ProposeWeight(fWeight);
320 if(aTrack.GetTrackStatus() != fAlive) { return theTotalResult; }
321
322 // Find cross section at end of step and check if <= 0
323 //
324 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
325 const G4Material* aMaterial = aTrack.GetMaterial();
326
327 // check only for charged particles
328 if(fXSType != fHadNoIntegral) {
331 theCrossSectionDataStore->ComputeCrossSection(aParticle,aMaterial);
332 //G4cout << "xs=" << xs << " xs0=" << theLastCrossSection
333 // << " " << aMaterial->GetName() << G4endl;
335 // No interaction
336 return theTotalResult;
337 }
338 }
339
340 const G4Element* anElement =
341 theCrossSectionDataStore->SampleZandA(aParticle,aMaterial,targetNucleus);
342
343 // Next check for illegal track status
344 //
345 if (aTrack.GetTrackStatus() != fAlive &&
346 aTrack.GetTrackStatus() != fSuspend) {
347 if (aTrack.GetTrackStatus() == fStopAndKill ||
351 ed << "G4HadronicProcess: track in unusable state - "
352 << aTrack.GetTrackStatus() << G4endl;
353 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
354 DumpState(aTrack,"PostStepDoIt",ed);
355 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
356 }
357 // No warning for fStopButAlive which is a legal status here
358 return theTotalResult;
359 }
360
361 // Initialize the hadronic projectile from the track
362 thePro.Initialise(aTrack);
363
364 theInteraction = ChooseHadronicInteraction(thePro, targetNucleus,
365 aMaterial, anElement);
366 if(nullptr == theInteraction) {
368 ed << "Target element "<<anElement->GetName()<<" Z= "
369 << targetNucleus.GetZ_asInt() << " A= "
370 << targetNucleus.GetA_asInt() << G4endl;
371 DumpState(aTrack,"ChooseHadronicInteraction",ed);
372 ed << " No HadronicInteraction found out" << G4endl;
373 G4Exception("G4HadronicProcess::PostStepDoIt", "had005",
374 FatalException, ed);
375 return theTotalResult;
376 }
377
378 G4HadFinalState* result = nullptr;
379 G4int reentryCount = 0;
380 /*
381 G4cout << "### " << aParticle->GetDefinition()->GetParticleName()
382 << " Ekin(MeV)= " << aParticle->GetKineticEnergy()
383 << " Z= " << targetNucleus.GetZ_asInt()
384 << " A= " << targetNucleus.GetA_asInt()
385 << " by " << theInteraction->GetModelName()
386 << G4endl;
387 */
388 do
389 {
390 try
391 {
392 // Call the interaction
393 result = theInteraction->ApplyYourself( thePro, targetNucleus);
394 ++reentryCount;
395 }
396 catch(G4HadronicException & aR)
397 {
399 aR.Report(ed);
400 ed << "Call for " << theInteraction->GetModelName() << G4endl;
401 ed << "Target element "<<anElement->GetName()<<" Z= "
402 << targetNucleus.GetZ_asInt()
403 << " A= " << targetNucleus.GetA_asInt() << G4endl;
404 DumpState(aTrack,"ApplyYourself",ed);
405 ed << " ApplyYourself failed" << G4endl;
406 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
407 ed);
408 }
409
410 // Check the result for catastrophic energy non-conservation
411 result = CheckResult(thePro, targetNucleus, result);
412
413 if(reentryCount>100) {
415 ed << "Call for " << theInteraction->GetModelName() << G4endl;
416 ed << "Target element "<<anElement->GetName()<<" Z= "
417 << targetNucleus.GetZ_asInt()
418 << " A= " << targetNucleus.GetA_asInt() << G4endl;
419 DumpState(aTrack,"ApplyYourself",ed);
420 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
421 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
422 ed);
423 }
424 }
425 while(!result); /* Loop checking, 30-Oct-2015, G.Folger */
426
427 // Check whether kaon0 or anti_kaon0 are present between the secondaries:
428 // if this is the case, transform them into either kaon0S or kaon0L,
429 // with equal, 50% probability, keeping their dynamical masses (and
430 // the other kinematical properties).
431 // When this happens - very rarely - a "JustWarning" exception is thrown.
432 // Because Fluka-Cern produces kaon0 and anti_kaon0, we reduce the number
433 // of warnings to max 1 per thread.
434 G4int nSec = (G4int)result->GetNumberOfSecondaries();
435 if ( nSec > 0 ) {
436 for ( G4int i = 0; i < nSec; ++i ) {
437 auto dynamicParticle = result->GetSecondary(i)->GetParticle();
438 auto part = dynamicParticle->GetParticleDefinition();
439 if ( part == G4KaonZero::Definition() ||
440 part == G4AntiKaonZero::Definition() ) {
441 G4ParticleDefinition* newPart;
442 if ( G4UniformRand() > 0.5 ) { newPart = G4KaonZeroShort::Definition(); }
443 else { newPart = G4KaonZeroLong::Definition(); }
444 dynamicParticle->SetDefinition( newPart );
445 }
446 }
447 }
448
449 result->SetTrafoToLab(thePro.GetTrafoToLab());
450 FillResult(result, aTrack);
451
452 if (epReportLevel != 0) {
453 CheckEnergyMomentumConservation(aTrack, targetNucleus);
454 }
455 //G4cout << "PostStepDoIt done nICelectrons= " << nICelectrons << G4endl;
456 return theTotalResult;
457}
@ FatalException
@ fKillTrackAndSecondaries
@ fSuspend
@ fPostponeToNextEvent
#define G4UniformRand()
Definition Randomize.hh:52
static G4AntiKaonZero * Definition()
const G4ParticleDefinition * GetParticleDefinition() const
const G4String & GetName() const
Definition G4Element.hh:115
void SetTrafoToLab(const G4LorentzRotation &aT)
void Report(std::ostream &aS) const
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4HadProjectile thePro
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
static G4KaonZeroLong * Definition()
static G4KaonZeroShort * Definition()
static G4KaonZero * Definition()
Definition G4KaonZero.cc:48
G4TrackStatus GetTrackStatus() const
G4double GetWeight() const
G4double theNumberOfInteractionLengthLeft

◆ PostStepGetPhysicalInteractionLength()

G4double G4HadronicProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Reimplemented in G4HadronStoppingProcess, G4MuNeutrinoNucleusProcess, G4NeutrinoElectronProcess, G4NeutronGeneralProcess, and G4TauNeutrinoNucleusProcess.

Definition at line 263 of file G4HadronicProcess.cc.

267{
269
270 const G4Material* mat = track.GetMaterial();
271 if(mat != currentMat) {
272 currentMat = mat;
274 matIdx = (G4int)track.GetMaterial()->GetIndex();
275 }
276 UpdateCrossSectionAndMFP(track.GetKineticEnergy());
277
278 // zero cross section
279 if(theLastCrossSection <= 0.0) {
282 return DBL_MAX;
283 }
284
285 // non-zero cross section
289 } else {
291 previousStepSize/currentInteractionLength;
294 }
297}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double G4Log(G4double x)
Definition G4Log.hh:169
std::size_t GetIndex() const
G4double currentInteractionLength
G4double theInitialNumberOfInteractionLength

◆ PreparePhysicsTable()

void G4HadronicProcess::PreparePhysicsTable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4HadronStoppingProcess, and G4NeutronGeneralProcess.

Definition at line 165 of file G4HadronicProcess.cc.

166{
167 if(nullptr == firstParticle) { firstParticle = &p; }
168 theProcessStore->RegisterParticle(this, &p);
169}

◆ ProcessDescription()

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

◆ RegisterMe()

void G4HadronicProcess::RegisterMe ( G4HadronicInteraction * a)

Definition at line 136 of file G4HadronicProcess.cc.

137{
138 if(nullptr == a) { return; }
139 theEnergyRangeManager.RegisterMe( a );
141}
void RegisterInteraction(G4HadronicProcess *, G4HadronicInteraction *)

Referenced by G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPiKBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4NeutronPHPBuilder::Build(), G4NeutronPHPBuilder::Build(), G4NeutronPHPBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4QGSPProtonBuilder::Build(), G4TritonPHPBuilder::Build(), G4HadProcesses::BuildNeutronInelasticAndCapture(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructProcess(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronElasticPhysicsHP::ConstructProcess(), G4HadronElasticPhysicsHPT::ConstructProcess(), G4HadronElasticPhysicsLEND::ConstructProcess(), G4HadronElasticPhysicsPHP::ConstructProcess(), G4HadronElasticPhysicsVI::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4HadronInelasticQBBC::ConstructProcess(), G4HadronInelasticQBBC_ABLA::ConstructProcess(), G4HadronPhysicsFTFQGSP_BERT::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4ThermalNeutrons::ConstructProcess(), G4URRNeutrons::ConstructProcess(), G4HadronicAbsorptionBertini::G4HadronicAbsorptionBertini(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4HadronicAbsorptionINCLXX::G4HadronicAbsorptionINCLXX(), G4MuonMinusCapture::G4MuonMinusCapture(), G4HadronPhysicsFTF_BIC::Neutron(), G4HadronPhysicsFTFP_BERT::Neutron(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), and G4HadronPhysicsShielding::Neutron().

◆ SetCrossSectionType()

void G4HadronicProcess::SetCrossSectionType ( G4HadXSType val)
inline

Definition at line 296 of file G4HadronicProcess.hh.

297{
298 fXSType = val;
299}

◆ SetEnergyMomentumCheckLevels()

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

Definition at line 317 of file G4HadronicProcess.hh.

319{
320 epCheckLevels.first = relativeLevel;
321 epCheckLevels.second = absoluteLevel;
322 levelsSetByProcess = true;
323}

◆ SetEpReportLevel()

void G4HadronicProcess::SetEpReportLevel ( G4int level)
inline

Definition at line 311 of file G4HadronicProcess.hh.

312{
313 epReportLevel = level;
314}

◆ SetIntegral()

void G4HadronicProcess::SetIntegral ( G4bool val)
inline

Definition at line 306 of file G4HadronicProcess.hh.

307{
308 useIntegralXS = val;
309}

◆ StartTracking()

void G4HadronicProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4NeutronGeneralProcess.

Definition at line 255 of file G4HadronicProcess.cc.

256{
257 currentMat = nullptr;
258 currentParticle = track->GetDefinition();
259 fDynParticle = track->GetDynamicParticle();
261}

◆ TwoPeaksXS()

std::vector< G4TwoPeaksHadXS * > * G4HadronicProcess::TwoPeaksXS ( ) const
inline

Definition at line 338 of file G4HadronicProcess.hh.

339{
340 return fXSpeaks;
341}

Member Data Documentation

◆ aScaleFactor

G4double G4HadronicProcess::aScaleFactor = 1.0
protected

◆ epReportLevel

◆ fWeight

G4double G4HadronicProcess::fWeight = 1.0
protected

Definition at line 220 of file G4HadronicProcess.hh.

Referenced by FillResult(), and PostStepDoIt().

◆ fXSType

◆ mfpKinEnergy

G4double G4HadronicProcess::mfpKinEnergy = DBL_MAX
protected

◆ theCrossSectionDataStore

◆ theLastCrossSection

G4double G4HadronicProcess::theLastCrossSection = 0.0
protected

◆ thePro

◆ theTotalResult


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