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

#include <G4VXRayModel.hh>

Inheritance diagram for G4VXRayModel:

Public Member Functions

 G4VXRayModel (const G4String &nam)
 G4VXRayModel (const G4VXRayModel &)
virtual ~G4VXRayModel ()
G4double Initialise (std::vector< const G4LogicalVolume * > *)
virtual void InitialiseModel ()
G4bool StepLimit (std::size_t idx, const G4Track &, G4double preStepBeta, G4double &limit)
virtual G4bool StepLimitForVolume (G4double &limit)
virtual void SampleXRays (std::vector< G4Track * > &out, const G4Step &)
virtual void ModelDescription (std::ostream &outFile) const
const G4StringGetName () const
G4int GetType () const
G4VXRayModeloperator= (const G4VXRayModel &right)=delete

Protected Attributes

std::vector< const G4LogicalVolume * > * pLogicalVolumes {nullptr}
G4LossTableManagerpEmManager {nullptr}
const G4LogicalVolumepCurrentLV {nullptr}
const G4TrackpCurrentTrack {nullptr}
G4double pBetaMin {1.0}
G4double pPreStepBeta {0.0}
G4double pMaxBetaChange {0.1}
G4int pMaxPhotons {100}
std::size_t pIndex {0}
std::size_t nVolumes {0}
G4int pTypeInt {0}
G4int pVerbose {1}
G4bool isMaster {false}
const G4String pName

Detailed Description

Definition at line 60 of file G4VXRayModel.hh.

Constructor & Destructor Documentation

◆ G4VXRayModel() [1/2]

G4VXRayModel::G4VXRayModel ( const G4String & nam)
explicit

Definition at line 51 of file G4VXRayModel.cc.

52 : pName(nam)
53{
54 Register();
55}
const G4String pName

Referenced by G4StandardCerenkovModel::G4StandardCerenkovModel(), G4VXRayModel(), and operator=().

◆ G4VXRayModel() [2/2]

G4VXRayModel::G4VXRayModel ( const G4VXRayModel & right)

Definition at line 59 of file G4VXRayModel.cc.

60 : pLogicalVolumes(nullptr),
61 pCurrentLV(nullptr),
62 pCurrentTrack(nullptr),
63 pBetaMin(1.0),
64 pPreStepBeta(0.0),
65 pMaxBetaChange(0.1),
66 pMaxPhotons(100),
67 pTypeInt(0),
68 pVerbose(1),
69 isMaster(false),
70 pName(right.pName)
71{
72 Register();
73}
std::vector< const G4LogicalVolume * > * pLogicalVolumes
const G4LogicalVolume * pCurrentLV
G4double pMaxBetaChange
G4double pBetaMin
const G4Track * pCurrentTrack
G4double pPreStepBeta

◆ ~G4VXRayModel()

G4VXRayModel::~G4VXRayModel ( )
virtual

Definition at line 77 of file G4VXRayModel.cc.

78{
79 pEmManager->DeRegister(this);
80}
G4LossTableManager * pEmManager

Member Function Documentation

◆ GetName()

const G4String & G4VXRayModel::GetName ( ) const
inline

Definition at line 90 of file G4VXRayModel.hh.

90{ return pName; };

Referenced by G4GeneralCerenkov::AddModelForVolume(), and G4LossTableManager::Register().

◆ GetType()

G4int G4VXRayModel::GetType ( ) const
inline

Definition at line 92 of file G4VXRayModel.hh.

92{ return pTypeInt; };

◆ Initialise()

G4double G4VXRayModel::Initialise ( std::vector< const G4LogicalVolume * > * ptr)

Definition at line 93 of file G4VXRayModel.cc.

94{
95 // definition of logical volumes
96 pLogicalVolumes = ptr;
97 if (nullptr != ptr) { nVolumes = (G4int)(ptr->size()); }
98 auto params = G4OpticalParameters::Instance();
99
100 // these parameters used by Cerenkov models
101 pMaxBetaChange = params->GetCerenkovMaxBetaChange();
102 pMaxPhotons = params->GetCerenkovMaxPhotonsPerStep();
103
104 // common initialisation
105 pVerbose = params->GetCerenkovVerboseLevel();
107
108 // needed for Cerenkov process
109 return pBetaMin;
110}
int G4int
Definition G4Types.hh:85
static G4OpticalParameters * Instance()
virtual void InitialiseModel()
std::size_t nVolumes

◆ InitialiseModel()

void G4VXRayModel::InitialiseModel ( )
virtual

Reimplemented in G4StandardCerenkovModel.

Definition at line 114 of file G4VXRayModel.cc.

115{}

Referenced by Initialise().

◆ ModelDescription()

void G4VXRayModel::ModelDescription ( std::ostream & outFile) const
virtual

Reimplemented in G4StandardCerenkovModel.

Definition at line 147 of file G4VXRayModel.cc.

148{
149 outFile << "The description for this model has not been written yet.\n";
150}

◆ operator=()

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

◆ SampleXRays()

void G4VXRayModel::SampleXRays ( std::vector< G4Track * > & out,
const G4Step &  )
virtual

Reimplemented in G4StandardCerenkovModel.

Definition at line 142 of file G4VXRayModel.cc.

143{}

◆ StepLimit()

G4bool G4VXRayModel::StepLimit ( std::size_t idx,
const G4Track & track,
G4double preStepBeta,
G4double & limit )

Definition at line 119 of file G4VXRayModel.cc.

121{
122 if (preStepBeta <= pBetaMin) { return false; }
123 pCurrentLV = nullptr;
124 if (idx < nVolumes) {
125 pIndex = idx;
126 pCurrentLV = (*pLogicalVolumes)[idx];
127 }
128 pCurrentTrack = &track;
129 pPreStepBeta = preStepBeta;
130 return StepLimitForVolume(limit);
131}
std::size_t pIndex
virtual G4bool StepLimitForVolume(G4double &limit)

◆ StepLimitForVolume()

G4bool G4VXRayModel::StepLimitForVolume ( G4double & limit)
virtual

Reimplemented in G4StandardCerenkovModel.

Definition at line 135 of file G4VXRayModel.cc.

136{
137 return false;
138}

Referenced by StepLimit().

Member Data Documentation

◆ isMaster

G4bool G4VXRayModel::isMaster {false}
protected

Definition at line 117 of file G4VXRayModel.hh.

117{false};

Referenced by G4VXRayModel().

◆ nVolumes

std::size_t G4VXRayModel::nVolumes {0}
protected

Definition at line 113 of file G4VXRayModel.hh.

113{0};

Referenced by Initialise(), G4StandardCerenkovModel::InitialiseModel(), and StepLimit().

◆ pBetaMin

G4double G4VXRayModel::pBetaMin {1.0}
protected

Definition at line 107 of file G4VXRayModel.hh.

107{1.0};

Referenced by G4VXRayModel(), Initialise(), G4StandardCerenkovModel::InitialiseModel(), and StepLimit().

◆ pCurrentLV

const G4LogicalVolume* G4VXRayModel::pCurrentLV {nullptr}
protected

Definition at line 105 of file G4VXRayModel.hh.

105{nullptr};

Referenced by G4VXRayModel(), and StepLimit().

◆ pCurrentTrack

const G4Track* G4VXRayModel::pCurrentTrack {nullptr}
protected

Definition at line 106 of file G4VXRayModel.hh.

106{nullptr};

Referenced by G4VXRayModel(), StepLimit(), and G4StandardCerenkovModel::StepLimitForVolume().

◆ pEmManager

G4LossTableManager* G4VXRayModel::pEmManager {nullptr}
protected

Definition at line 104 of file G4VXRayModel.hh.

104{nullptr};

Referenced by ~G4VXRayModel().

◆ pIndex

std::size_t G4VXRayModel::pIndex {0}
protected

◆ pLogicalVolumes

std::vector<const G4LogicalVolume*>* G4VXRayModel::pLogicalVolumes {nullptr}
protected

Definition at line 103 of file G4VXRayModel.hh.

103{nullptr};

Referenced by G4VXRayModel(), and Initialise().

◆ pMaxBetaChange

G4double G4VXRayModel::pMaxBetaChange {0.1}
protected

◆ pMaxPhotons

G4int G4VXRayModel::pMaxPhotons {100}
protected

Definition at line 111 of file G4VXRayModel.hh.

111{100};

Referenced by G4VXRayModel(), Initialise(), and G4StandardCerenkovModel::StepLimitForVolume().

◆ pName

const G4String G4VXRayModel::pName
protected

Definition at line 118 of file G4VXRayModel.hh.

Referenced by G4VXRayModel(), G4VXRayModel(), and GetName().

◆ pPreStepBeta

G4double G4VXRayModel::pPreStepBeta {0.0}
protected

◆ pTypeInt

G4int G4VXRayModel::pTypeInt {0}
protected

Definition at line 115 of file G4VXRayModel.hh.

115{0};

Referenced by G4VXRayModel(), and GetType().

◆ pVerbose

G4int G4VXRayModel::pVerbose {1}
protected

Definition at line 116 of file G4VXRayModel.hh.

116{1};

Referenced by G4VXRayModel(), and Initialise().


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