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

#include <G4PositronNuclearProcess.hh>

Inheritance diagram for G4PositronNuclearProcess:

Public Member Functions

 G4PositronNuclearProcess (const G4String &processName="positronNuclear")
virtual ~G4PositronNuclearProcess ()
virtual void ProcessDescription (std::ostream &outFile) const
Public Member Functions inherited from G4HadronInelasticProcess
 G4HadronInelasticProcess (const G4String &processName, const G4ParticleDefinition *particleDef=nullptr)
 ~G4HadronInelasticProcess () override
G4bool IsApplicable (const G4ParticleDefinition &) override
const G4ParticleDefinitionGetParticleDefinition () const
Public Member Functions inherited from G4HadronicProcess
 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 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 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 &)

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
Protected Member Functions inherited from G4HadronicProcess
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 inherited from G4HadronicProcess
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

Detailed Description

Definition at line 40 of file G4PositronNuclearProcess.hh.

Constructor & Destructor Documentation

◆ G4PositronNuclearProcess()

G4PositronNuclearProcess::G4PositronNuclearProcess ( const G4String & processName = "positronNuclear")

Definition at line 35 of file G4PositronNuclearProcess.cc.

38{
39 G4CrossSectionDataStore * theStore = GetCrossSectionDataStore();
40 theStore->AddDataSet(new G4ElectroNuclearCrossSection);
41}
void AddDataSet(G4VCrossSectionDataSet *)
G4HadronInelasticProcess(const G4String &processName, const G4ParticleDefinition *particleDef=nullptr)
G4CrossSectionDataStore * GetCrossSectionDataStore()
static G4Positron * Positron()
Definition G4Positron.cc:90

◆ ~G4PositronNuclearProcess()

G4PositronNuclearProcess::~G4PositronNuclearProcess ( )
virtual

Definition at line 44 of file G4PositronNuclearProcess.cc.

45{}

Member Function Documentation

◆ ProcessDescription()

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

Reimplemented from G4HadronicProcess.

Definition at line 48 of file G4PositronNuclearProcess.cc.

49{
50 outFile << "G4PositronNuclearProcess handles inelastic positron scattering\n"
51 << "from nuclei by invoking one hybrid electromagnetic-hadronic\n"
52 << "model and one hybrid electromagnetic-hadronic cross section set.\n";
53}

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