Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNASancheExcitationModel.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//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
29
30// Created by Z. Francis
31
32#ifndef G4DNASancheExcitationModel_h
33#define G4DNASancheExcitationModel_h 1
34
35#include <deque>
37
38#include "G4VEmModel.hh"
40#include "G4Electron.hh"
41#include "G4NistManager.hh"
42
44{
45public:
47 const G4String& nam = "DNASancheExcitationModel");
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 // Cross section
67
70
71 inline void ExtendLowEnergyLimit(G4double /*threshold*/);
72
73 inline void SetVerboseLevel(G4int verbose)
74 {
75 verboseLevel = verbose;
76 }
77
78 inline void SelectStationary(G4bool input);
79
80protected:
81
83
84private:
85
86 G4bool statCode;
87
88 // Water density table
89 const std::vector<G4double>* fpWaterDensity;
90
91 G4bool isInitialised;
92 G4int verboseLevel;
93
94 // Cross section
95
96 G4int RandomSelect(G4double energy);
97 G4int nLevels;
98 G4double VibrationEnergy(G4int level);
99 G4double Sum(G4double k);
100 G4double LinInterpolate(G4double e1,
101 G4double e2,
102 G4double e,
103 G4double xs1,
104 G4double xs2);
105
106 //
107// typedef std::map<double, std::map<double, double> > TriDimensionMap;
108// TriDimensionMap map1;
109 std::vector<G4double> tdummyVec;
110 std::vector<std::vector<G4double>> fEnergyLevelXS;
111 std::vector<G4double> fEnergyTotalXS;
112
113 //
116
117};
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
122{
123 if(threshold < 2 * CLHEP::eV)
124 G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not "
125 "validated below 2 eV !",
126 "", JustWarning, "");
127
128 SetLowEnergyLimit(threshold);
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
134{
135 statCode = input;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
140#endif
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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 void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4double PartialCrossSection(G4double energy, G4int level)
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753