Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAPTBAugerModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Authors: S. Meylan and C. Villagrasa (IRSN, France)
27// Models come from
28// M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
29//
30//
31// -------------------------------------------------------------------
32//
33// Geant4 Header G4DNAPTBAugerModel
34//
35// -------------------------------------------------------------------
36//
37// Class description:
38// Implementation of atomic deexcitation
39//
40// -------------------------------------------------------------------
41
42#ifndef G4DNAPTBAugerModel_h
43#define G4DNAPTBAugerModel_h 1
44
46#include "G4AtomicShell.hh"
47#include "globals.hh"
48#include "G4DynamicParticle.hh"
49#include <vector>
50
53class G4EmCorrections;
54class G4Material;
55
56/*!
57 * \brief The G4DNAPTBAugerModel class
58 * Implement the PTB Auger model
59 */
61{
62public:
63
64 /*!
65 * \brief G4DNAPTBAugerModel
66 * Constructor
67 * \param modelName
68 */
69 G4DNAPTBAugerModel(const G4String &modelName);
70
71 /*!
72 * \brief ~G4DNAPTBAugerModel
73 * Destructor
74 */
75 virtual ~G4DNAPTBAugerModel();
76
77 G4DNAPTBAugerModel(G4DNAPTBAugerModel &) = delete; // prevent copy-construction
78 G4DNAPTBAugerModel & operator=(const G4DNAPTBAugerModel &right) = delete; // prevent assignement
79
80 /*!
81 * \brief Initialise
82 * Set the verbose value
83 */
84 virtual void Initialise();
85
86 /*!
87 * \brief SetCutForAugerElectrons
88 * Set the cut for the auger electrons production
89 * \param cut
90 */
92
93 /*!
94 * \brief ComputeAugerEffect
95 * Main method to be called by the ionisation model.
96 * \param fvect
97 * \param materialNameIni
98 * \param bindingEnergy
99 */
100 void ComputeAugerEffect(std::vector<G4DynamicParticle *> *fvect, const G4String& materialNameIni, G4double bindingEnergy);
101
102private:
103
104 const G4String modelName; ///< name of the auger model
105
106 G4int verboseLevel;
107 G4double minElectronEnergy;
108
109 /*!
110 * \brief GenerateAugerWithRandomDirection
111 * Generates the auger particle
112 * \param fvect
113 * \param kineticEnergy
114 */
115 void GenerateAugerWithRandomDirection(std::vector<G4DynamicParticle*>* fvect, G4double kineticEnergy);
116
117 /*!
118 * \brief CalculAugerEnergyFor
119 * \param atomId
120 * \return the auger particle energy
121 */
122 G4double CalculAugerEnergyFor(G4int atomId);
123
124 /*!
125 * \brief DetermineIonisedAtom
126 * \param atomId
127 * \param materialName
128 * \param bindingEnergy
129 * \return the id of the chosen ionised atom
130 */
131 G4int DetermineIonisedAtom(G4int atomId, const G4String &materialName, G4double bindingEnergy);
132
133
134};
135
136#endif
137
138
139
140
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
The G4DNAPTBAugerModel class Implement the PTB Auger model.
G4DNAPTBAugerModel(G4DNAPTBAugerModel &)=delete
virtual ~G4DNAPTBAugerModel()
~G4DNAPTBAugerModel Destructor
void ComputeAugerEffect(std::vector< G4DynamicParticle * > *fvect, const G4String &materialNameIni, G4double bindingEnergy)
ComputeAugerEffect Main method to be called by the ionisation model.
void SetCutForAugerElectrons(G4double cut)
SetCutForAugerElectrons Set the cut for the auger electrons production.
G4DNAPTBAugerModel & operator=(const G4DNAPTBAugerModel &right)=delete
virtual void Initialise()
Initialise Set the verbose value.
G4DNAPTBAugerModel(const G4String &modelName)
G4DNAPTBAugerModel Constructor.