Geant4-11
G4DNADiracRMatrixExcitationModel.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// Created on 2016/04/08
27//
28// Authors: D. Sakata, S. Incerti
29//
30// This class perform electric excitation for electron transportation,
31// based on Dirac B-Spline R-Matrix Model and scaled experimental data.
32// See following reference paper
33// Phys.Rev.A77,062711(2008) and Phys.Rev.A78,042713(2008)
34
35#ifndef G4DNADiracRMatrixExcitationModel_h
36#define G4DNADiracRMatrixExcitationModel_h 1
37
38#include "G4VEmModel.hh"
42
44#include "G4Electron.hh"
45#include "G4Proton.hh"
46#include "G4NistManager.hh"
47
49
51{
52
53public:
54
56 const G4String& nam = "DNADiracRMatrixExcitationModel");
57
59
60 virtual void Initialise(const G4ParticleDefinition*,
61 const G4DataVector& = *(new G4DataVector()));
62
64 const G4ParticleDefinition* p,
65 G4double ekin,
66 G4double emin,
68
71 G4double kineticEnergy);
72
74 G4int level,
76 G4double kineticEnergy);
77
78 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
80 const G4DynamicParticle*,
81 G4double tmin,
82 G4double maxEnergy);
83
84 inline void SelectStationary(G4bool input);
85
86protected:
87
89
90private:
91
92 const G4double paramFuncTCS_5dto6s1[3]={-3e-50 , 9.46358e-16, 1.4237 }; // y = [0]+[1]/pow(x-[2],2)
93 const G4double paramFuncTCS_5dto6s2[3]={-3e-50 , 4.24498e-15, -0.674543}; // y = [0]+[1]/pow(x-[2],2)
94 const G4double paramFuncTCS_6sto6p1[3]={ 1.50018e-26, 2.459e-15 ,-40.8088 }; // y = [0]+[1]*log(x-[2])/(x-[2])
95 const G4double paramFuncTCS_6sto6p2[3]={ 1.26684e-25, 3.97221e-15,-55.6954 }; // y = [0]+[1]*log(x-[2])/(x-[2])
96 const G4int ShellEnumAu [4]={19 , 20 ,21 , 21 };
97 // 5d3/2 ,6s1/2 ,6s1/2 //from EADL
98 const G4double BindingEnergyAu [4]={12.16 ,10.46 , 8.3 , 8.3 };
99 // [eV] 5d3/2 ,6s1/2 ,6s1/2 //from EADL
100 const G4double ExcitationEnergyAu[4]={ 2.66 , 1.14 , 4.63 , 5.11};
101 // [eV] 5dto6s1,6sto6p1,6sto6p2
102
106
110
113 const std::vector<G4double>* fpMaterialDensity=nullptr;
116
119 G4double kineticEnergy);
120
121
123 =(const G4DNADiracRMatrixExcitationModel &right);
125
126};
127
129{
130 statCode = input;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
134
135#endif
static const G4double emax
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const std::vector< G4double > * fpMaterialDensity
G4DNADiracRMatrixExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNADiracRMatrixExcitationModel")
virtual G4double GetExtendedTotalCrossSection(const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
virtual G4double GetExtendedPartialCrossSection(const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
G4int RandomSelect(const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
G4DNADiracRMatrixExcitationModel(const G4DNADiracRMatrixExcitationModel &)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
string material
Definition: eplot.py:19