00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052 #ifndef G4eSingleCoulombScatteringModel_h
00053 #define G4eSingleCoulombScatteringModel_h 1
00054
00055 #include "G4VEmModel.hh"
00056 #include "globals.hh"
00057 #include "G4NistManager.hh"
00058 #include "G4ScreeningMottCrossSection.hh"
00059
00060 #include <vector>
00061 using namespace std;
00062
00063 class G4ParticleChangeForGamma;
00064 class G4ParticleDefinition;
00065
00066 class G4eSingleCoulombScatteringModel : public G4VEmModel
00067 {
00068
00069 public:
00070
00071 G4eSingleCoulombScatteringModel(const G4String& nam = "eSingleCoulombScattering");
00072
00073 virtual ~G4eSingleCoulombScatteringModel();
00074
00075 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00076
00077 virtual G4double ComputeCrossSectionPerAtom(
00078 const G4ParticleDefinition*,
00079 G4double kinEnergy,
00080 G4double Z,
00081 G4double A,
00082 G4double cut,
00083 G4double emax);
00084
00085
00086
00087
00088
00089
00090 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00091 const G4MaterialCutsCouple*,
00092 const G4DynamicParticle*,
00093 G4double tmin,
00094 G4double maxEnergy);
00095
00096 inline void SetRecoilThreshold(G4double eth);
00097
00098
00099
00100 protected:
00101
00102
00103 inline void DefineMaterial(const G4MaterialCutsCouple*);
00104 inline void SetupParticle(const G4ParticleDefinition*);
00105
00106
00107
00108 private:
00109
00110
00111
00112 G4eSingleCoulombScatteringModel & operator=(const G4eSingleCoulombScatteringModel &right);
00113 G4eSingleCoulombScatteringModel(const G4eSingleCoulombScatteringModel&);
00114
00115
00116 protected:
00117
00118
00119 G4ParticleTable* theParticleTable;
00120 G4ParticleChangeForGamma* fParticleChange;
00121 G4NistManager* fNistManager;
00122 G4ScreeningMottCrossSection* Mottcross;
00123
00124 const std::vector<G4double>* pCuts;
00125 const G4MaterialCutsCouple* currentCouple;
00126 const G4Material* currentMaterial;
00127 const G4Element* currentElement;
00128 G4int currentMaterialIndex;
00129
00130
00131 G4double cosThetaMin;
00132 G4double recoilThreshold;
00133
00134
00135
00136 const G4ParticleDefinition* particle;
00137 G4double mass;
00138 G4double lowEnergyLimit;
00139
00140
00141 private:
00142
00143 G4bool isInitialised;
00144
00145 };
00146
00147
00148 void G4eSingleCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
00149 {
00150 if(cup != currentCouple) {
00151 currentCouple = cup;
00152 currentMaterial = cup->GetMaterial();
00153 currentMaterialIndex = currentCouple->GetIndex();
00154
00155 }
00156 }
00157
00158
00159
00160 inline
00161 void G4eSingleCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
00162 {
00163 if(p != particle) {
00164 particle = p;
00165 mass = particle->GetPDGMass();
00166 Mottcross->SetupParticle(p);
00167 }
00168 }
00169
00170
00171 inline void G4eSingleCoulombScatteringModel::SetRecoilThreshold(G4double eth)
00172 {
00173 recoilThreshold = eth;
00174 }
00175
00176
00177
00178 #endif