Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNARPWBAIonisationModel.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// Reference:
27// A.D. Dominguez-Munoz, M.I. Gallardo, M.C. Bordage,
28// Z. Francis, S. Incerti, M.A. Cortes-Giraldo,
29// Radiat. Phys. Chem. 199 (2022) 110363.
30//
31// Class authors:
32// A.D. Dominguez-Munoz
33// M.A. Cortes-Giraldo (miancortes -at- us.es)
34//
35// Class creation: 2022-03-03
36//
37//
38
39#ifndef G4DNARPWBAIonisationModel_h
40#define G4DNARPWBAIonisationModel_h 1
41#include "G4VEmModel.hh"
45#include "G4Electron.hh"
46#include "G4Proton.hh"
51#include "G4NistManager.hh"
52
54{
55 public:
57 const G4String& nam = "DNARPWBAIonisationModel");
60 delete;
63 const G4DataVector& = *(new G4DataVector())) override;
64
66 const G4ParticleDefinition* p, G4double ekin,
67 G4double emin, G4double emax) override;
68 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
70 G4double tmin, G4double maxEnergy) override;
73 G4double /*kineticEnergy*/) override;
74 G4double DifferentialCrossSection(const G4double& k, const G4double& energyTransfer,
75 const G4int& shell);
76 G4double TransferedEnergy(G4double incomingParticleEnergy, G4int shell,
77 const G4double& random);
78 inline void SelectFasterComputation(const G4bool& input);
79 inline void SelectStationary(const G4bool& input);
80 inline void SelectSPScaling(const G4bool& input);
81
82 protected:
84
85 private:
86 using MapData = std::map<G4String, std::unique_ptr<G4DNACrossSectionDataSet>,
87 std::less<G4String>>;
88 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
89 using VecMap = std::map<G4double, std::vector<G4double>>;
90 // methods
91 G4double RandomizeEjectedElectronEnergy(const G4double&, const G4int& );
92 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(const G4double& incomingParticleEnergy,
93 const G4int& shell);
94 G4double Interpolate(const G4double& e1, const G4double& e2, const G4double& e, const G4double& xs1,
95 const G4double& xs2);
96 G4double QuadInterpolator(const G4double& e11, const G4double& e12, const G4double& e21,
97 const G4double& e22, const G4double& x11, const G4double& x12,
98 const G4double& x21, const G4double& x22, const G4double& t1,
99 const G4double& t2, const G4double& t, const G4double& e);
100 //shoud change to array ?
101 G4int RandomSelect(G4double energy);
102 void InitialiseForProton(const G4ParticleDefinition*);
103 G4bool InEnergyLimit(const G4double&);
104
105 //members
106 G4bool fasterCode = false;
107 G4bool statCode = false;
108 G4bool spScaling = true;
109 // Water density table
110 const std::vector<G4double>* fpMolWaterDensity = nullptr;
111 // Deexcitation manager to produce fluo photons and e-
112 G4VAtomDeexcitation* fAtomDeexcitation = nullptr;
113 G4double lowEnergyLimit = 0;
114 G4double highEnergyLimit = 0;
115 G4bool isInitialised = false;
116 G4int verboseLevel = 0;
117 // Cross section
118
119 std::unique_ptr<G4DNACrossSectionDataSet> fpTotalCrossSection;
120 // Final state
121 G4DNAWaterIonisationStructure waterStructure;
122 TriDimensionMap eDiffCrossSectionData[6];
123 TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
124 TriDimensionMap pDiffCrossSectionData[6];
125 TriDimensionMap pNrjTransfData[6]; // for cumulated dcs
126 std::vector<G4double> eTdummyVec;
127 std::vector<G4double> pTdummyVec;
128 VecMap eVecm;
129 VecMap pVecm;
130 VecMap eProbaShellMap[6]; // for cumulated dcs
131 VecMap pProbaShellMap[6]; // for cumulated dcs
133};
134
136 const G4bool& input)
137{
138 fasterCode = input;
139}
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
143{
144 statCode = input;
145}
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147
149{
150 spScaling = input;
151}
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153
154#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DNARPWBAIonisationModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNARPWBAIonisationModel")
void SelectFasterComputation(const G4bool &input)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNARPWBAIonisationModel(const G4DNARPWBAIonisationModel &)=delete
void SelectStationary(const G4bool &input)
void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
void SelectSPScaling(const G4bool &input)
G4DNARPWBAIonisationModel & operator=(const G4DNARPWBAIonisationModel &right)=delete
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double TransferedEnergy(G4double incomingParticleEnergy, G4int shell, const G4double &random)
G4double DifferentialCrossSection(const G4double &k, const G4double &energyTransfer, const G4int &shell)
G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double) override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85