Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAIonElasticModel.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// Author: H. N. Tran (Ton Duc Thang University)
27// p, H, He, He+ and He++ models are assumed identical
28// NIMB 343, 132-137 (2015)
29//
30// The Geant4-DNA web site is available at http://geant4-dna.org
31//
32
33#ifndef G4DNAIonElasticModel_h
34#define G4DNAIonElasticModel_h 1
35
36#include <map>
38#include "G4VEmModel.hh"
39#include "G4Proton.hh"
44#include "G4NistManager.hh"
45
47{
48
49public:
50
52 const G4String& nam ="DNAIonElasticModel");
53
54
55 ~G4DNAIonElasticModel () override;
56
58 operator= (const G4DNAIonElasticModel &right) = delete;
60
61 void
62 Initialise (const G4ParticleDefinition* particuleDefinition,
63 const G4DataVector&) override;
64
66 CrossSectionPerVolume (const G4Material* material,
67 const G4ParticleDefinition* p, G4double ekin,
68 G4double emin, G4double emax) override;
69
70 void
71 SampleSecondaries (std::vector<G4DynamicParticle*>*,
73 G4double tmin, G4double maxEnergy) override;
74
75 void
77
80 {
81 return killBelowEnergy;
82 }
83
84 inline void SelectStationary(G4bool input);
85
86protected:
87
89
90private:
91
92 G4bool statCode;
93
94 // Water density table
95 const std::vector<G4double>* fpMolWaterDensity;
96
97 G4double killBelowEnergy;
98 G4double lowEnergyLimit;
99 G4double highEnergyLimit;
100 G4bool isInitialised{false};
101 G4int verboseLevel;
102
103 G4double fParticle_Mass;
104
105 // Cross section
106 G4DNACrossSectionDataSet* fpTableData;
107
108 // Final state
109
111 Theta (G4ParticleDefinition * aParticleDefinition, G4double k,
112 G4double integrDiff);
113
115 LinLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1,
116 G4double xs2);
117
119 LogLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1,
120 G4double xs2);
121
123 QuadInterpolator (G4double e11, G4double e12, G4double e21, G4double e22,
124 G4double x11, G4double x12, G4double x21, G4double x22,
125 G4double t1, G4double t2, G4double t, G4double e);
126
128 LinLinInterpolate (G4double e1, G4double e2, G4double e, G4double xs1,
129 G4double xs2);
130
131 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
132 TriDimensionMap fDiffCrossSectionData;
133
134 std::vector<G4double> eTdummyVec;
135
136 using VecMap = std::map<G4double, std::vector<G4double>>;
137 VecMap eVecm;
138
140 RandomizeThetaCM (G4double k, G4ParticleDefinition * aParticleDefinition);
141
142};
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145
147{
148 statCode = input;
149}
150
151#endif
152
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetKillBelowThreshold(G4double threshold)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4DNAIonElasticModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAIonElasticModel")
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
void Initialise(const G4ParticleDefinition *particuleDefinition, const G4DataVector &) override
G4DNAIonElasticModel(const G4DNAIonElasticModel &)=delete
G4DNAIonElasticModel & operator=(const G4DNAIonElasticModel &right)=delete
void SelectStationary(G4bool input)
G4ParticleChangeForGamma * fParticleChangeForGamma