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

#include <G4NIELCalculator.hh>

Public Member Functions

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

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:67
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
void SetNIELCalculator(G4NIELCalculator *)
const G4String & GetName() const

◆ ~G4NIELCalculator()

G4NIELCalculator::~G4NIELCalculator ( )
default

◆ G4NIELCalculator() [2/2]

G4NIELCalculator::G4NIELCalculator ( const G4NIELCalculator & )
delete

Member Function Documentation

◆ AddEmModel()

void G4NIELCalculator::AddEmModel ( G4VEmModel * mod)

Definition at line 72 of file G4NIELCalculator.cc.

73{
74 if(mod && mod != fModel) {
75 fModel = mod;
76 if(fVerbose > 0) {
77 G4cout << "G4NIELCalculator: new model <" << fModel->GetName()
78 << "> is added" << G4endl;
79 }
80 }
81}

◆ ComputeNIEL()

G4double G4NIELCalculator::ComputeNIEL ( const G4Step * step)

Definition at line 90 of file G4NIELCalculator.cc.

91{
92 G4double niel = 0.0;
94 if(fModel && T2 > 0.) {
95 const G4Track* track = step->GetTrack();
96 const G4ParticleDefinition* part = track->GetParticleDefinition();
97 G4double length = step->GetStepLength();
98
99 if(length > 0.0 && part->GetPDGMass() > 100*CLHEP::MeV) {
100
101 // primary
103 G4double T = 0.5*(T1 + T2);
104 const G4MaterialCutsCouple* couple =
106 niel = length*fModel->ComputeDEDXPerVolume(couple->GetMaterial(),part,T);
107 niel = std::min(niel, T1);
108 }
109 }
110 return niel;
111}
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)

◆ Initialise()

void G4NIELCalculator::Initialise ( )

Definition at line 85 of file G4NIELCalculator.cc.

86{}

Referenced by G4LossTableManager::BuildPhysicsTable().

◆ operator=()

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

◆ RecoilEnergy()

G4double G4NIELCalculator::RecoilEnergy ( const G4Step * step)

Definition at line 115 of file G4NIELCalculator.cc.

116{
117 G4double erec = 0.0;
118 const std::vector<const G4Track*>* sec = step->GetSecondaryInCurrentStep();
119
120 if(sec) {
121 for(auto track : *sec) {
122 const G4ParticleDefinition* part = track->GetParticleDefinition();
123 if(part->IsGeneralIon()) {
124 erec += track->GetKineticEnergy();
125 }
126 }
127 }
128 return erec;
129}
G4bool IsGeneralIon() const
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition G4Step.cc:202

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