source: trunk/source/processes/electromagnetic/lowenergy/include/G4DNAChampionElasticModel.hh @ 1350

Last change on this file since 1350 was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 4.9 KB
Line 
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,v 1.4 2010/11/11 22:32:22 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29
30#ifndef G4DNAChampionElasticModel_h
31#define G4DNAChampionElasticModel_h 1
32
33#include <map>
34#include "G4DNACrossSectionDataSet.hh"
35#include "G4VEmModel.hh"
36#include "G4Electron.hh"
37#include "G4ParticleChangeForGamma.hh"
38#include "G4LogLogInterpolation.hh"
39#include "G4ProductionCutsTable.hh"
40
41class G4DNAChampionElasticModel : public G4VEmModel
42{
43
44public:
45
46  G4DNAChampionElasticModel(const G4ParticleDefinition* p = 0, 
47                            const G4String& nam = "DNAChampionElasticModel");
48
49  virtual ~G4DNAChampionElasticModel();
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  inline void SetKillBelowThreshold (G4double threshold);               
66  G4double GetKillBelowThreshold () { return killBelowEnergy; }         
67
68protected:
69
70  G4ParticleChangeForGamma* fParticleChangeForGamma;
71
72private:
73
74  G4double killBelowEnergy; 
75  G4double lowEnergyLimit; 
76  G4double highEnergyLimit; 
77  G4bool isInitialised;
78  G4int verboseLevel;
79 
80  // Cross section
81 
82  typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
83  MapFile tableFile;
84
85  typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
86  MapData tableData;
87 
88  // Final state
89
90  //G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double theta);
91
92  G4double Theta(G4ParticleDefinition * aParticleDefinition, G4double k, G4double integrDiff);
93 
94  G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
95
96  G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
97   
98  G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
99   
100  G4double QuadInterpolator(G4double e11, 
101                            G4double e12, 
102                            G4double e21, 
103                            G4double e22, 
104                            G4double x11,
105                            G4double x12, 
106                            G4double x21, 
107                            G4double x22, 
108                            G4double t1, 
109                            G4double t2, 
110                            G4double t, 
111                            G4double e);
112
113  typedef std::map<double, std::map<double, double> > TriDimensionMap;
114
115  TriDimensionMap eDiffCrossSectionData;
116  std::vector<double> eTdummyVec;
117
118  typedef std::map<double, std::vector<double> > VecMap;
119  VecMap eVecm;
120   
121  G4double RandomizeCosTheta(G4double k);
122   
123  //
124   
125  G4DNAChampionElasticModel & operator=(const  G4DNAChampionElasticModel &right);
126  G4DNAChampionElasticModel(const  G4DNAChampionElasticModel&);
127
128};
129
130
131inline void G4DNAChampionElasticModel::SetKillBelowThreshold (G4double threshold) 
132{ 
133    killBelowEnergy = threshold; 
134   
135    if (threshold < 1*eV)
136     G4Exception ("*** WARNING : the G4DNAChampionElasticModel class is not validated below 1 eV !","",JustWarning,"") ;
137   
138    if (threshold < 0.025*eV) threshold = 0.025*eV;
139             
140}               
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
143
144#endif
Note: See TracBrowser for help on using the repository browser.