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

The G4DNAPTBAugerModel class Implement the PTB Auger model. More...

#include <G4DNAPTBAugerModel.hh>

Public Member Functions

void ComputeAugerEffect (std::vector< G4DynamicParticle * > *fvect, const G4String &materialNameIni, G4double bindingEnergy)
 ComputeAugerEffect Main method to be called by the ionisation model. More...
 
 G4DNAPTBAugerModel (const G4String &modelName)
 G4DNAPTBAugerModel Constructor. More...
 
virtual void Initialise ()
 Initialise Set the verbose value. More...
 
void SetCutForAugerElectrons (G4double cut)
 SetCutForAugerElectrons Set the cut for the auger electrons production. More...
 
virtual ~G4DNAPTBAugerModel ()
 ~G4DNAPTBAugerModel Destructor More...
 

Private Member Functions

G4double CalculAugerEnergyFor (G4int atomId)
 CalculAugerEnergyFor. More...
 
G4int DetermineIonisedAtom (G4int atomId, const G4String &materialName, G4double bindingEnergy)
 DetermineIonisedAtom. More...
 
 G4DNAPTBAugerModel (G4DNAPTBAugerModel &)
 
void GenerateAugerWithRandomDirection (std::vector< G4DynamicParticle * > *fvect, G4double kineticEnergy)
 GenerateAugerWithRandomDirection Generates the auger particle. More...
 
G4DNAPTBAugerModeloperator= (const G4DNAPTBAugerModel &right)
 

Private Attributes

G4double minElectronEnergy
 
const G4String modelName
 name of the auger model More...
 
G4int verboseLevel
 

Detailed Description

The G4DNAPTBAugerModel class Implement the PTB Auger model.

Definition at line 60 of file G4DNAPTBAugerModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBAugerModel() [1/2]

G4DNAPTBAugerModel::G4DNAPTBAugerModel ( const G4String modelName)

G4DNAPTBAugerModel Constructor.

Parameters
modelName

Definition at line 41 of file G4DNAPTBAugerModel.cc.

41 : modelName(modelAugerName)
42{
43 // To inform the user that the Auger model is enabled
44 G4cout << modelName <<" is constructed" << G4endl;
45}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String modelName
name of the auger model

References G4cout, G4endl, and modelName.

◆ ~G4DNAPTBAugerModel()

G4DNAPTBAugerModel::~G4DNAPTBAugerModel ( )
virtual

~G4DNAPTBAugerModel Destructor

Definition at line 47 of file G4DNAPTBAugerModel.cc.

48{
49 if( verboseLevel>0 ) G4cout << modelName <<" is deleted" << G4endl;
50}

References G4cout, G4endl, modelName, and verboseLevel.

◆ G4DNAPTBAugerModel() [2/2]

G4DNAPTBAugerModel::G4DNAPTBAugerModel ( G4DNAPTBAugerModel )
private

Member Function Documentation

◆ CalculAugerEnergyFor()

G4double G4DNAPTBAugerModel::CalculAugerEnergyFor ( G4int  atomId)
private

CalculAugerEnergyFor.

Parameters
atomId
Returns
the auger particle energy

Definition at line 136 of file G4DNAPTBAugerModel.cc.

137{
138 G4double kineticEnergy;
139
140 if(atomId==2) // oxygen
141 {
142 kineticEnergy = 495*eV;
143 }
144 else
145 {
146 G4double f1, f2, f3, g1, g2, Y;
147
148 Y = G4UniformRand();
149
150 if(atomId == 1){ // carbon
151 f1 = -7.331e-2;
152 f2 = -3.306e-5;
153 f3 = 2.433e0;
154 g1 = 4.838e-1;
155 g2 = 3.886e0;
156 }
157 else if(atomId == 4){ // nitrogen
158 f1 = -7.518e-2;
159 f2 = 1.178e-4;
160 f3 = 2.600e0;
161 g1 = 4.639e-1;
162 g2 = 3.770e0;
163 }
164 else// if(atomId == 3) // carbon_TMP
165 {
166 f1 = -5.700e-2;
167 f2 = 1.200e-4;
168 f3 = 2.425e0;
169 g1 = 5.200e-1;
170 g2 = 2.560e0;
171 }
172
173 kineticEnergy = pow(10, f1*pow( abs( log10(Y) ) , g1) + f2*pow( abs( log10(Y) ) , g2) + f3 )*eV;
174 }
175
176 return kineticEnergy;
177}
G4double Y(G4double density)
static constexpr double eV
Definition: G4SIunits.hh:201
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52

References eV, G4UniformRand, and Y().

Referenced by ComputeAugerEffect().

◆ ComputeAugerEffect()

void G4DNAPTBAugerModel::ComputeAugerEffect ( std::vector< G4DynamicParticle * > *  fvect,
const G4String materialNameIni,
G4double  bindingEnergy 
)

ComputeAugerEffect Main method to be called by the ionisation model.

Parameters
fvect
materialNameIni
bindingEnergy

Definition at line 65 of file G4DNAPTBAugerModel.cc.

66{
67 // Rename material if modified NIST material
68 // This is needed when material is obtained from G4MaterialCutsCouple
69 G4String materialName = materialNameIni;
70 if(materialName.find("_MODIFIED")){
71 materialName = materialName.substr(0,materialName.size()-9);
72 }
73
74 // check if there is a k-shell ionisation and find the ionised atom
75 G4int atomId(0);
76
77 atomId = DetermineIonisedAtom(atomId, materialName, bindingEnergy);
78
79 if(atomId!=0)
80 {
81 G4double kineticEnergy = CalculAugerEnergyFor(atomId);
82
83 if(kineticEnergy<0)
84 {
85 G4cerr<<"**************************"<<G4endl;
86 G4cerr<<"FatalError. Auger kineticEnergy: "<<kineticEnergy<<G4endl;
87 exit(EXIT_FAILURE);
88 }
89
90 if(atomId==1 || atomId==2 || atomId==3)
91 {
92 GenerateAugerWithRandomDirection(fvect, kineticEnergy);
93 }
94 else if(atomId==4)
95 {
96 GenerateAugerWithRandomDirection(fvect, kineticEnergy);
97 GenerateAugerWithRandomDirection(fvect, kineticEnergy);
98 }
99 }
100}
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
void GenerateAugerWithRandomDirection(std::vector< G4DynamicParticle * > *fvect, G4double kineticEnergy)
GenerateAugerWithRandomDirection Generates the auger particle.
G4int DetermineIonisedAtom(G4int atomId, const G4String &materialName, G4double bindingEnergy)
DetermineIonisedAtom.
G4double CalculAugerEnergyFor(G4int atomId)
CalculAugerEnergyFor.
G4double bindingEnergy(G4int A, G4int Z)
def exit()
Definition: g4zmq.py:99

References G4InuclSpecialFunctions::bindingEnergy(), CalculAugerEnergyFor(), DetermineIonisedAtom(), g4zmq::exit(), G4cerr, G4endl, and GenerateAugerWithRandomDirection().

Referenced by G4DNAPTBIonisationModel::SampleSecondaries().

◆ DetermineIonisedAtom()

G4int G4DNAPTBAugerModel::DetermineIonisedAtom ( G4int  atomId,
const G4String materialName,
G4double  bindingEnergy 
)
private

DetermineIonisedAtom.

Parameters
atomId
materialName
bindingEnergy
Returns
the id of the chosen ionised atom

Definition at line 104 of file G4DNAPTBAugerModel.cc.

105{
106 if(materialName=="THF" || materialName=="backbone_THF"){
107 if(bindingEnergy==305.07){
108 atomId=1; //"carbon";
109 }
110 else if(bindingEnergy==557.94){
111 atomId=2; //"oxygen";
112 }
113 }
114 else if(materialName=="PY" || materialName=="PU"
115 || materialName=="cytosine_PY" || materialName=="thymine_PY"
116 || materialName=="adenine_PU" || materialName=="guanine_PU"
117 )
118 {
119 if(bindingEnergy==307.52){
120 atomId=1; //"carbon";
121 }
122 else if(bindingEnergy==423.44){
123 atomId=4; //"nitrogen";
124 }
125 }
126 else if(materialName=="TMP"|| materialName=="backbone_TMP"){
127 if(bindingEnergy==209.59 || bindingEnergy==152.4)
128 atomId=3; //"carbonTMP";
129 }
130
131 return atomId;
132}

References G4InuclSpecialFunctions::bindingEnergy().

Referenced by ComputeAugerEffect().

◆ GenerateAugerWithRandomDirection()

void G4DNAPTBAugerModel::GenerateAugerWithRandomDirection ( std::vector< G4DynamicParticle * > *  fvect,
G4double  kineticEnergy 
)
private

GenerateAugerWithRandomDirection Generates the auger particle.

Parameters
fvect
kineticEnergy

Definition at line 188 of file G4DNAPTBAugerModel.cc.

189{
190 // Isotropic angular distribution for the outcoming e-
191 G4double newcosTh = 1.-2.*G4UniformRand();
192 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
193 G4double newPhi = twopi*G4UniformRand();
194
195 G4double xDir = newsinTh*std::sin(newPhi);
196 G4double yDir = newsinTh*std::cos(newPhi);
197 G4double zDir = newcosTh;
198
199 G4ThreeVector ElectronDirection(xDir,yDir,zDir);
200
201 // generation of new particle
202 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(), ElectronDirection, kineticEnergy) ;
203 fvect->push_back(dp);
204}
static constexpr double twopi
Definition: G4SIunits.hh:56
static G4Electron * Electron()
Definition: G4Electron.cc:93

References G4Electron::Electron(), G4UniformRand, and twopi.

Referenced by ComputeAugerEffect().

◆ Initialise()

void G4DNAPTBAugerModel::Initialise ( )
virtual

Initialise Set the verbose value.

Definition at line 52 of file G4DNAPTBAugerModel.cc.

53{
54 verboseLevel = 0;
55
56 if( verboseLevel>0 )
57 {
58 G4cout << "PTB Auger model is initialised " << G4endl;
59 }
60
61}

References G4cout, G4endl, and verboseLevel.

Referenced by G4DNAPTBIonisationModel::Initialise().

◆ operator=()

G4DNAPTBAugerModel & G4DNAPTBAugerModel::operator= ( const G4DNAPTBAugerModel right)
private

◆ SetCutForAugerElectrons()

void G4DNAPTBAugerModel::SetCutForAugerElectrons ( G4double  cut)

SetCutForAugerElectrons Set the cut for the auger electrons production.

Parameters
cut

Definition at line 181 of file G4DNAPTBAugerModel.cc.

182{
183 minElectronEnergy = cut;
184}

References minElectronEnergy.

Field Documentation

◆ minElectronEnergy

G4double G4DNAPTBAugerModel::minElectronEnergy
private

Definition at line 105 of file G4DNAPTBAugerModel.hh.

Referenced by SetCutForAugerElectrons().

◆ modelName

const G4String G4DNAPTBAugerModel::modelName
private

name of the auger model

Definition at line 102 of file G4DNAPTBAugerModel.hh.

Referenced by G4DNAPTBAugerModel(), and ~G4DNAPTBAugerModel().

◆ verboseLevel

G4int G4DNAPTBAugerModel::verboseLevel
private

Definition at line 104 of file G4DNAPTBAugerModel.hh.

Referenced by Initialise(), and ~G4DNAPTBAugerModel().


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