Geant4-11
G4DNAELSEPAElasticModel.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// Created on 2016/01/18
27//
28// Authors: D. Sakata, S. Incerti
29//
30// Based on a recent release of the ELSEPA code
31// developed and provided kindly by F. Salvat et al.
32// See
33// Computer Physics Communications, 165(2), 157-190. (2005)
34// http://dx.doi.org/10.1016/j.cpc.2004.09.006
35//
36
37#ifndef G4DNAELSEPAElasticModel_h
38#define G4DNAELSEPAElasticModel_h 1
39
40#include <map>
42#include "G4VEmModel.hh"
43#include "G4Electron.hh"
47#include "G4NistManager.hh"
48
50{
51
52public:
53
55 const G4String& nam = "DNAELSEPAElasticModel");
56
58
59 virtual void Initialise(
60 const G4ParticleDefinition* particle, const G4DataVector&);
61
63 const G4ParticleDefinition* particle,
64 G4double ekin,
65 G4double emin,
67
68 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
70 const G4DynamicParticle*,
71 G4double tmin,
72 G4double maxEnergy);
73
75 {fhighEnergyLimit = input; SetHighEnergyLimit(input);};
76
77 void SetKillBelowThreshold (G4double threshold);
78
80
81protected:
82
84
85private:
86
87 const std::vector<G4double>* fpMolDensity=nullptr;
88 std::vector <G4double> kIntersectionEnergySR;
89
93
96
99
100 G4double Theta(G4int Z, G4ParticleDefinition * aParticleDefinition,
101 G4double k,
102 G4double integrDiff);
103
105 G4double e2,
106 G4double e,
107 G4double xs1,
108 G4double xs2);
109
111 G4double e2,
112 G4double e,
113 G4double xs1,
114 G4double xs2);
115
117 G4double e2,
118 G4double e,
119 G4double xs1,
120 G4double xs2);
121
123 G4double e2,
124 G4double e,
125 G4double xs1,
126 G4double xs2);
127
129 G4double e12,
130 G4double e21,
131 G4double e22,
132 G4double x11,
133 G4double x12,
134 G4double x21,
135 G4double x22,
136 G4double t1,
137 G4double t2,
138 G4double t,
139 G4double e);
140
142
143 typedef std::map<G4int,std::map<G4double,std::map<G4double,G4double>>>
146
147 std::map <G4int, std::vector<G4double> > eEdummyVecZ;
148
149 typedef std::map<G4double, std::vector<G4double> > VecMap;
152
153 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
156
157 std::vector<G4double> eEdummyVec_Au;
158 std::vector<G4double> eEdummyVec_H2O;
159
162
163};
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167#endif
static const G4double e1[44]
static const G4double e2[44]
static const G4double emax
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
std::map< G4int, std::map< G4double, std::map< G4double, G4double > > > TriDimensionMapZ
const std::vector< G4double > * fpMolDensity
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *particle, G4double ekin, G4double emin, G4double emax)
std::map< G4double, std::map< G4double, G4double > > TriDimensionMap
std::vector< G4double > eEdummyVec_Au
G4double LogLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
std::map< G4int, std::vector< G4double > > eEdummyVecZ
G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4double Theta(G4int Z, G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
G4DNAELSEPAElasticModel(const G4DNAELSEPAElasticModel &)
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
std::map< G4double, std::vector< G4double > > VecMap
G4DNAELSEPAElasticModel(const G4ParticleDefinition *particle=0, const G4String &nam="DNAELSEPAElasticModel")
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
std::vector< G4double > eEdummyVec_H2O
void SetMaximumEnergy(G4double input)
G4DNAELSEPAElasticModel & operator=(const G4DNAELSEPAElasticModel &right)
std::vector< G4double > kIntersectionEnergySR
void SetKillBelowThreshold(G4double threshold)
G4double RandomizeCosTheta(G4int Z, G4double k)
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNACrossSectionDataSet * fpData_H2O
G4DNACrossSectionDataSet * fpData_Au
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
string material
Definition: eplot.py:19