Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAChampionElasticModel.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//
27
28#ifndef G4DNAChampionElasticModel_h
29#define G4DNAChampionElasticModel_h 1
30
31#include <map>
33#include "G4VEmModel.hh"
34#include "G4Electron.hh"
38#include "G4NistManager.hh"
39
41{
42
43public:
44
46 const G4String& nam = "DNAChampionElasticModel");
47
49
50 virtual void Initialise(const G4ParticleDefinition*,
51 const G4DataVector&);
52
53 virtual G4double CrossSectionPerVolume(const G4Material* material,
54 const G4ParticleDefinition* p,
55 G4double ekin,
56 G4double emin,
57 G4double emax);
58
59 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
61 const G4DynamicParticle*,
62 G4double tmin,
63 G4double maxEnergy);
64
65 void SetKillBelowThreshold(G4double threshold);
66
68 {
70 errMsg << "The method G4DNAChampionElasticModel::"
71 "GetKillBelowThreshold is deprecated";
72
73 G4Exception("G4DNAChampionElasticModel::GetKillBelowThreshold",
74 "deprecated",
76 errMsg);
77 return 0.;
78 }
79
80private:
81 // Cross section
82 typedef std::map<G4double, std::vector<G4double> > VecMap;
83 VecMap eVecm;
84 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
85 TriDimensionMap eDiffCrossSectionData;
86 std::vector<G4double> eTdummyVec;
87
88 // Water density table
89 const std::vector<G4double>* fpMolWaterDensity;
90
91 // Cross section
93
94protected:
96
97private:
98 G4int verboseLevel;
99 G4bool isInitialised;
100
101 // Final state
102
103 //G4double DifferentialCrossSection(G4ParticleDefinition* aParticle,
104 // G4double k, G4double theta);
105
106 G4double Theta(//G4ParticleDefinition * aParticleDefinition,
107 G4double k,
108 G4double integrDiff);
109
110 G4double LinLinInterpolate(G4double e1,
111 G4double e2,
112 G4double e,
113 G4double xs1,
114 G4double xs2);
115
116 G4double LinLogInterpolate(G4double e1,
117 G4double e2,
118 G4double e,
119 G4double xs1,
120 G4double xs2);
121
122 G4double LogLogInterpolate(G4double e1,
123 G4double e2,
124 G4double e,
125 G4double xs1,
126 G4double xs2);
127
128 G4double QuadInterpolator(G4double e11,
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
141 G4double RandomizeCosTheta(G4double k);
142
143 //
144
147};
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150
151#endif
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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
void SetKillBelowThreshold(G4double threshold)
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)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChangeForGamma