Geant4  10.02
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 // $Id: G4DNAChampionElasticModel.hh 87137 2014-11-25 09:12:48Z gcosmo $
27 //
28 
29 #ifndef G4DNAChampionElasticModel_h
30 #define G4DNAChampionElasticModel_h 1
31 
32 #include <map>
34 #include "G4VEmModel.hh"
35 #include "G4Electron.hh"
37 #include "G4LogLogInterpolation.hh"
38 #include "G4ProductionCutsTable.hh"
39 #include "G4NistManager.hh"
40 
42 {
43 
44 public:
45 
47  const G4String& nam = "DNAChampionElasticModel");
48 
50 
51  virtual void Initialise(const G4ParticleDefinition*, 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*>*,
60  const G4MaterialCutsCouple*,
61  const G4DynamicParticle*,
62  G4double tmin,
63  G4double maxEnergy);
64 
65  void SetKillBelowThreshold(G4double threshold);
66 
68  {
69  return killBelowEnergy;
70  }
71 
72 protected:
73 
75 
76 private:
77  // Water density table
78  const std::vector<G4double>* fpMolWaterDensity;
79 
85 
86  // Cross section
87 
88  typedef std::map<G4String, G4String, std::less<G4String> > MapFile;
89  MapFile tableFile;
90 
91  typedef std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> > MapData;
92  MapData tableData;
93 
94  // Final state
95 
96  //G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double theta);
97 
98  G4double Theta(G4ParticleDefinition * aParticleDefinition,
99  G4double k,
100  G4double integrDiff);
101 
103  G4double e2,
104  G4double e,
105  G4double xs1,
106  G4double xs2);
107 
109  G4double e2,
110  G4double e,
111  G4double xs1,
112  G4double xs2);
113 
115  G4double e2,
116  G4double e,
117  G4double xs1,
118  G4double xs2);
119 
121  G4double e12,
122  G4double e21,
123  G4double e22,
124  G4double x11,
125  G4double x12,
126  G4double x21,
127  G4double x22,
128  G4double t1,
129  G4double t2,
130  G4double t,
131  G4double e);
132 
133  typedef std::map<double, std::map<double, double> > TriDimensionMap;
134 
135  TriDimensionMap eDiffCrossSectionData;
136  std::vector<double> eTdummyVec;
137 
138  typedef std::map<double, std::vector<double> > VecMap;
139  VecMap eVecm;
140 
142 
143  //
144 
147 
148 };
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
152 #endif
G4DNAChampionElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAChampionElasticModel")
void SetKillBelowThreshold(G4double threshold)
std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > MapData
static const G4double e2
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
int G4int
Definition: G4Types.hh:78
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
std::map< double, std::map< double, double > > TriDimensionMap
bool G4bool
Definition: G4Types.hh:79
G4double Theta(G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
static const G4double e1
static const G4double emax
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)
G4DNAChampionElasticModel & operator=(const G4DNAChampionElasticModel &right)
G4ParticleChangeForGamma * fParticleChangeForGamma
double G4double
Definition: G4Types.hh:76
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
std::map< double, std::vector< double > > VecMap
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
std::map< G4String, G4String, std::less< G4String > > MapFile
const std::vector< G4double > * fpMolWaterDensity