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

#include <G4OpAbsorption.hh>

+ Inheritance diagram for G4OpAbsorption:

Public Member Functions

 G4OpAbsorption (const G4String &processName="OpAbsorption", G4ProcessType type=fOptical)
 
virtual ~G4OpAbsorption ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void Initialise ()
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 
- 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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (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
 
virtual void ProcessDescription (std::ostream &outfile) 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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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 52 of file G4OpAbsorption.hh.

Constructor & Destructor Documentation

◆ G4OpAbsorption()

G4OpAbsorption::G4OpAbsorption ( const G4String processName = "OpAbsorption",
G4ProcessType  type = fOptical 
)
explicit

Definition at line 56 of file G4OpAbsorption.cc.

57 : G4VDiscreteProcess(processName, type)
58{
59 Initialise();
60 if(verboseLevel > 0)
61 {
62 G4cout << GetProcessName() << " is created " << G4endl;
63 }
65}
@ fOpAbsorption
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void Initialise()
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4OpAbsorption()

G4OpAbsorption::~G4OpAbsorption ( )
virtual

Definition at line 68 of file G4OpAbsorption.cc.

68{}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4OpAbsorption::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 102 of file G4OpAbsorption.cc.

104{
105 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
108 G4double attLength = DBL_MAX;
109
110 if(MPT)
111 {
113 if(attVector)
114 {
115 attLength =
116 attVector->Value(aParticle->GetTotalMomentum(), idx_absorption);
117 }
118 }
119
120 return attLength;
121}
double G4double
Definition: G4Types.hh:83
G4double GetTotalMomentum() const
G4MaterialPropertyVector * GetProperty(const char *key, G4bool warning=false)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:254
G4double Value(G4double theEnergy, std::size_t &lastidx) const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition: templates.hh:62

◆ Initialise()

void G4OpAbsorption::Initialise ( )
virtual

Definition at line 77 of file G4OpAbsorption.cc.

78{
79 SetVerboseLevel(G4OpticalParameters::Instance()->GetAbsorptionVerboseLevel());
80}
static G4OpticalParameters * Instance()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412

Referenced by G4OpAbsorption(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4OpAbsorption::IsApplicable ( const G4ParticleDefinition aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 85 of file G4OpAbsorption.hh.

87{
88 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
89}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 83 of file G4OpAbsorption.cc.

85{
87
88 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
89 G4double thePhotonMomentum = aParticle->GetTotalMomentum();
90
93
94 if(verboseLevel > 1)
95 {
96 G4cout << "\n** OpAbsorption: Photon absorbed! **" << G4endl;
97 }
98 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
99}
@ fStopAndKill
virtual void Initialize(const G4Track &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327

◆ PreparePhysicsTable()

void G4OpAbsorption::PreparePhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 71 of file G4OpAbsorption.cc.

72{
73 Initialise();
74}

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