Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// $Id: G4DNAELSEPAElasticModel.hh 97497 2016-06-03 11:41:57Z matkara $
27//
28
29#ifndef G4DNAELSEPAElasticModel_h
30#define G4DNAELSEPAElasticModel_h 1
31
32#include <map>
34#include "G4VEmModel.hh"
35#include "G4Electron.hh"
39#include "G4NistManager.hh"
40
42{
43
44public:
45
47 const G4String& nam = "DNAELSEPAElasticModel");
48
50
51 virtual void Initialise(const G4ParticleDefinition*,
52 const G4DataVector&);
53
54 virtual G4double CrossSectionPerVolume(const G4Material* material,
55 const G4ParticleDefinition* p,
56 G4double ekin,
57 G4double emin,
58 G4double emax);
59
60 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
62 const G4DynamicParticle*,
63 G4double tmin,
64 G4double maxEnergy);
65
66 void SetKillBelowThreshold(G4double threshold);
67
69 {
71 errMsg << "The method G4DNAELSEPAElasticModel::"
72 "GetKillBelowThreshold is deprecated";
73
74 G4Exception("G4DNAELSEPAElasticModel::GetKillBelowThreshold",
75 "deprecated",
77 errMsg);
78 return 0.;
79 }
80
81private:
82 // Cross section
83 typedef std::map<double, std::vector<double> > VecMap;
84 VecMap eVecm;
85 typedef std::map<double, std::map<double, double> > TriDimensionMap;
86 TriDimensionMap eDiffCrossSectionData;
87 std::vector<double> eTdummyVec;
88
89 // Water density table
90 const std::vector<G4double>* fpMolWaterDensity;
91
92 // Cross section
94
95protected:
97
98private:
99 G4int verboseLevel;
100 G4bool isInitialised;
101
102 // Final state
103
104 //G4double DifferentialCrossSection(G4ParticleDefinition* aParticle,
105 // G4double k, G4double theta);
106
107 G4double Theta(//G4ParticleDefinition * aParticleDefinition,
108 G4double k,
109 G4double integrDiff);
110
111 G4double LinLinInterpolate(G4double e1,
112 G4double e2,
113 G4double e,
114 G4double xs1,
115 G4double xs2);
116
117 G4double LinLogInterpolate(G4double e1,
118 G4double e2,
119 G4double e,
120 G4double xs1,
121 G4double xs2);
122
123 G4double LogLogInterpolate(G4double e1,
124 G4double e2,
125 G4double e,
126 G4double xs1,
127 G4double xs2);
128
129 G4double QuadInterpolator(G4double e11,
130 G4double e12,
131 G4double e21,
132 G4double e22,
133 G4double x11,
134 G4double x12,
135 G4double x21,
136 G4double x22,
137 G4double t1,
138 G4double t2,
139 G4double t,
140 G4double e);
141
142 G4double RandomizeCosTheta(G4double k);
143
144 //
145
146 G4DNAELSEPAElasticModel & operator=(const G4DNAELSEPAElasticModel &right);
148};
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151
152#endif
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void SetKillBelowThreshold(G4double threshold)
G4ParticleChangeForGamma * fParticleChangeForGamma