Geant4-11
Public Member Functions | Private Attributes
G4NIELCalculator Class Reference

#include <G4NIELCalculator.hh>

Public Member Functions

void AddEmModel (G4VEmModel *)
 
G4double ComputeNIEL (const G4Step *)
 
 G4NIELCalculator (const G4NIELCalculator &)=delete
 
 G4NIELCalculator (G4VEmModel *, G4int verb)
 
void Initialise ()
 
G4NIELCalculatoroperator= (const G4NIELCalculator &right)=delete
 
G4double RecoilEnergy (const G4Step *)
 
 ~G4NIELCalculator ()
 

Private Attributes

G4VEmModelfModel
 
G4int fVerbose
 

Detailed Description

Definition at line 60 of file G4NIELCalculator.hh.

Constructor & Destructor Documentation

◆ G4NIELCalculator() [1/2]

G4NIELCalculator::G4NIELCalculator ( G4VEmModel mod,
G4int  verb 
)
explicit

Definition at line 60 of file G4NIELCalculator.cc.

61 : fModel(mod), fVerbose(verb)
62{
64 if(fVerbose > 0) {
65 G4cout << "G4NIELCalculator: is created with the model <"
66 << fModel->GetName() << ">" << G4endl;
67 }
68}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
void SetNIELCalculator(G4NIELCalculator *)
G4VEmModel * fModel
const G4String & GetName() const
Definition: G4VEmModel.hh:837

References fModel, fVerbose, G4cout, G4endl, G4VEmModel::GetName(), G4LossTableManager::Instance(), and G4LossTableManager::SetNIELCalculator().

◆ ~G4NIELCalculator()

G4NIELCalculator::~G4NIELCalculator ( )

Definition at line 72 of file G4NIELCalculator.cc.

73{}

◆ G4NIELCalculator() [2/2]

G4NIELCalculator::G4NIELCalculator ( const G4NIELCalculator )
delete

Member Function Documentation

◆ AddEmModel()

void G4NIELCalculator::AddEmModel ( G4VEmModel mod)

Definition at line 77 of file G4NIELCalculator.cc.

78{
79 if(mod && mod != fModel) {
80 fModel = mod;
81 if(fVerbose > 0) {
82 G4cout << "G4NIELCalculator: new model <" << fModel->GetName()
83 << "> is added" << G4endl;
84 }
85 }
86}

References fModel, fVerbose, G4cout, G4endl, and G4VEmModel::GetName().

◆ ComputeNIEL()

G4double G4NIELCalculator::ComputeNIEL ( const G4Step step)

Definition at line 95 of file G4NIELCalculator.cc.

96{
97 G4double niel = 0.0;
99 if(fModel && T2 > 0.) {
100 const G4Track* track = step->GetTrack();
101 const G4ParticleDefinition* part = track->GetParticleDefinition();
102 G4double length = step->GetStepLength();
103
104 if(length > 0.0 && part->GetPDGMass() > 100*CLHEP::MeV) {
105
106 // primary
108 G4double T = 0.5*(T1 + T2);
109 const G4MaterialCutsCouple* couple =
111 niel = length*fModel->ComputeDEDXPerVolume(couple->GetMaterial(),part,T);
112 niel = std::min(niel, T1);
113 }
114 }
115 return niel;
116}
double G4double
Definition: G4Types.hh:83
const G4Material * GetMaterial() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:228
static constexpr double MeV
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4VEmModel::ComputeDEDXPerVolume(), fModel, G4StepPoint::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4StepPoint::GetMaterialCutsCouple(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetPDGMass(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4Step::GetStepLength(), G4Step::GetTrack(), CLHEP::MeV, and G4INCL::Math::min().

◆ Initialise()

void G4NIELCalculator::Initialise ( )

Definition at line 90 of file G4NIELCalculator.cc.

91{}

Referenced by G4LossTableManager::BuildPhysicsTable().

◆ operator=()

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

◆ RecoilEnergy()

G4double G4NIELCalculator::RecoilEnergy ( const G4Step step)

Definition at line 120 of file G4NIELCalculator.cc.

121{
122 G4double erec = 0.0;
123 const std::vector<const G4Track*>* sec = step->GetSecondaryInCurrentStep();
124
125 if(sec) {
126 for(auto track : *sec) {
127 const G4ParticleDefinition* part = track->GetParticleDefinition();
128 if(part->IsGeneralIon()) {
129 erec += track->GetKineticEnergy();
130 }
131 }
132 }
133 return erec;
134}
G4bool IsGeneralIon() const
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition: G4Step.cc:209

References G4Step::GetSecondaryInCurrentStep(), and G4ParticleDefinition::IsGeneralIon().

Field Documentation

◆ fModel

G4VEmModel* G4NIELCalculator::fModel
private

Definition at line 86 of file G4NIELCalculator.hh.

Referenced by AddEmModel(), ComputeNIEL(), and G4NIELCalculator().

◆ fVerbose

G4int G4NIELCalculator::fVerbose
private

Definition at line 87 of file G4NIELCalculator.hh.

Referenced by AddEmModel(), and G4NIELCalculator().


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