Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
78 /*!
79 * \brief Initialise
80 * Set the verbose value
81 */
82 virtual void Initialise();
83
84 /*!
85 * \brief SetCutForAugerElectrons
86 * Set the cut for the auger electrons production
87 * \param cut
88 */
90
91 /*!
92 * \brief ComputeAugerEffect
93 * Main method to be called by the ionisation model.
94 * \param fvect
95 * \param materialNameIni
96 * \param bindingEnergy
97 */
98 void ComputeAugerEffect(std::vector<G4DynamicParticle *> *fvect, const G4String& materialNameIni, G4double bindingEnergy);
99
100private:
101
102 const G4String modelName; ///< name of the auger model
103
104 G4int verboseLevel;
105 G4double minElectronEnergy;
106
107 /*!
108 * \brief GenerateAugerWithRandomDirection
109 * Generates the auger particle
110 * \param fvect
111 * \param kineticEnergy
112 */
113 void GenerateAugerWithRandomDirection(std::vector<G4DynamicParticle*>* fvect, G4double kineticEnergy);
114
115 /*!
116 * \brief CalculAugerEnergyFor
117 * \param atomId
118 * \return the auger particle energy
119 */
120 G4double CalculAugerEnergyFor(G4int atomId);
121
122 /*!
123 * \brief DetermineIonisedAtom
124 * \param atomId
125 * \param materialName
126 * \param bindingEnergy
127 * \return the id of the chosen ionised atom
128 */
129 G4int DetermineIonisedAtom(G4int atomId, const G4String &materialName, G4double bindingEnergy);
130
131 // copy constructor and hide assignment operator
132 G4DNAPTBAugerModel(G4DNAPTBAugerModel &); // prevent copy-construction
133 G4DNAPTBAugerModel & operator=(const G4DNAPTBAugerModel &right); // prevent assignement
134
135};
136
137#endif
138
139
140
141
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
The G4DNAPTBAugerModel class Implement the PTB Auger model.
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.
virtual void Initialise()
Initialise Set the verbose value.