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

#include <G4DNAPositronium.hh>

+ Inheritance diagram for G4DNAPositronium:

Public Member Functions

 G4DNAPositronium (const G4String &processName="G4DNAPositronium")
 
 ~G4DNAPositronium () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void InitialiseProcess (const G4ParticleDefinition *) override
 
virtual void PrintInfo ()
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
 ~G4VEmProcess () override
 
void ProcessDescription (std::ostream &outFile) const override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void StartTracking (G4Track *) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
virtual G4VEmProcessGetEmProcess (const G4String &name)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
 
G4double GetCrossSection (const G4double kinEnergy, const G4MaterialCutsCouple *couple) override
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
void SetLambdaBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
void SetLambdaTable (G4PhysicsTable *)
 
void SetLambdaTablePrim (G4PhysicsTable *)
 
std::vector< G4double > * EnergyOfCrossSectionMax () const
 
void SetEnergyOfCrossSectionMax (std::vector< G4double > *)
 
G4CrossSectionType CrossSectionType () const
 
void SetCrossSectionType (G4CrossSectionType val)
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t idxCouple) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4int NumberOfModels () const
 
G4VEmModelEmModel (size_t index=0) const
 
const G4VEmModelGetCurrentModel () const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4ElementGetCurrentElement () const
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
G4double CrossSectionBiasingFactor () const
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void SetEmMasterProcess (const G4VEmProcess *)
 
void SetBuildTableFlag (G4bool val)
 
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
 
G4bool UseBaseMaterial () const
 
void BuildLambdaTable ()
 
void StreamInfo (std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
 
 G4VEmProcess (G4VEmProcess &)=delete
 
G4VEmProcessoperator= (const G4VEmProcess &right)=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 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)
 
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 G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
 
G4VEmModelSelectModel (G4double kinEnergy, size_t)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
void DefineMaterial (const G4MaterialCutsCouple *couple)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
const G4MaterialCutsCoupleMaterialCutsCouple () const
 
G4bool ApplyCuts () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
G4int DensityIndex (G4int idx) const
 
G4double DensityFactor (G4int idx) const
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VEmProcess
const G4MaterialCutsCouplecurrentCouple = nullptr
 
const G4MaterialcurrentMaterial = nullptr
 
G4EmBiasingManagerbiasManager = nullptr
 
std::vector< G4double > * theEnergyOfCrossSectionMax = nullptr
 
G4double mfpKinEnergy = DBL_MAX
 
G4double preStepKinEnergy = 0.0
 
G4double preStepLambda = 0.0
 
G4int mainSecondaries = 1
 
G4int secID = _EM
 
G4int fluoID = _Fluorescence
 
G4int augerID = _AugerElectron
 
G4int biasID = _EM
 
G4int tripletID = _TripletElectron
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
size_t coupleIdxLambda = 0
 
size_t idxLambda = 0
 
G4bool isTheMaster = true
 
G4bool baseMat = false
 
std::vector< G4DynamicParticle * > secParticles
 
G4ParticleChangeForGamma fParticleChange
 
- 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 31 of file G4DNAPositronium.hh.

Constructor & Destructor Documentation

◆ G4DNAPositronium()

G4DNAPositronium::G4DNAPositronium ( const G4String & processName = "G4DNAPositronium")

Definition at line 31 of file G4DNAPositronium.cc.

32 : G4VEmProcess( processName )
33{
34} // constructor
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)

◆ ~G4DNAPositronium()

G4DNAPositronium::~G4DNAPositronium ( )
overridedefault

Member Function Documentation

◆ InitialiseProcess()

void G4DNAPositronium::InitialiseProcess ( const G4ParticleDefinition * )
overridevirtual

Implements G4VEmProcess.

Definition at line 48 of file G4DNAPositronium.cc.

49{
50 if(!isInitialised)
51 {
52 isInitialised = true;
53 SetBuildTableFlag(false);
54 if(EmModel() == nullptr) SetEmModel(new G4LEPTSPositroniumModel);
55 EmModel()->SetLowEnergyLimit(0.1*CLHEP::eV);
56 EmModel()->SetHighEnergyLimit(15.*CLHEP::MeV);
57 AddEmModel(1, EmModel());
58
59 }
60}
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
G4VEmModel * EmModel(size_t index=0) const
void SetBuildTableFlag(G4bool val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)

◆ IsApplicable()

G4bool G4DNAPositronium::IsApplicable ( const G4ParticleDefinition & particleDef)
overridevirtual

Implements G4VEmProcess.

Definition at line 42 of file G4DNAPositronium.cc.

42 {
43 return( &particleDef == G4Positron::Positron() );
44}
static G4Positron * Positron()
Definition G4Positron.cc:90

◆ PrintInfo()

void G4DNAPositronium::PrintInfo ( )
virtual

Definition at line 63 of file G4DNAPositronium.cc.

64{
65 G4cout
66 << " Total cross sections computed from " << EmModel()->GetName() << " model"
67 << G4endl;
68}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const

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