Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAPTBExcitationModel.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#ifndef G4DNAPTBExcitationModel_h
32#define G4DNAPTBExcitationModel_h 1
33
34#include "G4VDNAModel.hh"
37
40#include "G4Electron.hh"
41#include "G4Proton.hh"
42#include "G4NistManager.hh"
43
45
48
49/*!
50 * \brief The G4DNAPTBExcitationModel class
51 * This class implements the PTB excitation model.
52 */
54{
55
56public:
57
58 /*!
59 * \brief G4DNAPTBExcitationModel
60 * Constructor
61 * \param applyToMaterial
62 * \param p
63 * \param nam
64 */
65 G4DNAPTBExcitationModel(const G4String &applyToMaterial = "all", const G4ParticleDefinition* p = 0,
66 const G4String& nam = "DNAPTBExcitationModel");
67
68 /*!
69 * \brief ~G4DNAPTBExcitationModel
70 * Destructor
71 */
73
74 /*!
75 * \brief Initialise
76 * Set the materials for which the model can be used and defined the energy limits
77 */
78 virtual void Initialise(const G4ParticleDefinition* particle, const G4DataVector& = *(new G4DataVector()), G4ParticleChangeForGamma* fpChangeForGamme=nullptr);
79
80 /*!
81 * \brief CrossSectionPerVolume
82 * Retrieve the cross section corresponding to the current material, particle and energy
83 * \param material
84 * \param materialName
85 * \param p
86 * \param ekin
87 * \param emin
88 * \param emax
89 * \return the cross section value
90 */
91 virtual G4double CrossSectionPerVolume(const G4Material* material,
92 const G4String& materialName,
93 const G4ParticleDefinition* p,
94 G4double ekin,
95 G4double emin,
96 G4double emax);
97
98 /*!
99 * \brief SampleSecondaries
100 * If the model is selected for the ModelInterface then the SampleSecondaries method will be called.
101 * The method sets the incident particle characteristics after the ModelInterface.
102 * \param materialName
103 * \param particleChangeForGamma
104 * \param tmin
105 * \param tmax
106 */
107 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109 const G4String& materialName,
110 const G4DynamicParticle*,
111 G4ParticleChangeForGamma *particleChangeForGamma,
112 G4double tmin,
113 G4double tmax);
114
115protected:
116
117private:
118
119 G4int verboseLevel; ///< verbose level
120
121 G4DNAWaterExcitationStructure waterStructure;
122
123 G4DNAPTBExcitationStructure ptbExcitationStructure;
124 G4DNAPTBIonisationStructure ptbIonisationStructure;
125
126 typedef std::map<G4String,G4double,std::less<G4String> > MapMeanEnergy;
127 MapMeanEnergy tableMeanEnergyPTB; ///< map: [materialName]=energyValue
128
129 // copy constructor and hide assignment operator
130 G4DNAPTBExcitationModel(const G4DNAPTBExcitationModel&); // prevent copy-construction
131 G4DNAPTBExcitationModel & operator=(const G4DNAPTBExcitationModel &right); // prevent assignement
132};
133
134#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
The G4DNAPTBExcitationModel class This class implements the PTB excitation model.
virtual ~G4DNAPTBExcitationModel()
~G4DNAPTBExcitationModel Destructor
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume Retrieve the cross section corresponding to the current material,...
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &= *(new G4DataVector()), G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
Initialise Set the materials for which the model can be used and defined the energy limits.
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
SampleSecondaries If the model is selected for the ModelInterface then the SampleSecondaries method w...
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50