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
00053
00054
00055
00056
00057 #ifndef G4IonCoulombScatteringModel_h
00058 #define G4IonCoulombScatteringModel_h 1
00059
00060 #include "G4VEmModel.hh"
00061 #include "globals.hh"
00062 #include "G4NistManager.hh"
00063 #include "G4IonCoulombCrossSection.hh"
00064
00065 #include <vector>
00066 using namespace std;
00067
00068 class G4ParticleChangeForGamma;
00069 class G4ParticleDefinition;
00070
00071 class G4IonCoulombScatteringModel : public G4VEmModel
00072 {
00073
00074 public:
00075
00076 G4IonCoulombScatteringModel(const G4String& nam = "IonCoulombScattering");
00077
00078 virtual ~G4IonCoulombScatteringModel();
00079
00080 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00081
00082 virtual G4double ComputeCrossSectionPerAtom(
00083 const G4ParticleDefinition*,
00084 G4double kinEnergy,
00085 G4double Z,
00086 G4double A,
00087 G4double cut,
00088 G4double emax);
00089
00090 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00091 const G4MaterialCutsCouple*,
00092 const G4DynamicParticle*,
00093 G4double tmin,
00094 G4double maxEnergy);
00095
00096
00097
00098 inline void SetRecoilThreshold(G4double eth);
00099 void SetHeavyIonCorr(G4int b) {heavycorr=b; };
00100 G4int GetHeavyIonCorr() {return heavycorr; };
00101
00102
00103
00104
00105 protected:
00106
00107
00108 inline void DefineMaterial(const G4MaterialCutsCouple*);
00109
00110 inline void SetupParticle(const G4ParticleDefinition*);
00111
00112
00113
00114 private:
00115
00116
00117
00118 G4IonCoulombScatteringModel & operator=(const G4IonCoulombScatteringModel &right);
00119 G4IonCoulombScatteringModel(const G4IonCoulombScatteringModel&);
00120
00121
00122 protected:
00123
00124
00125 G4ParticleTable* theParticleTable;
00126 G4ParticleChangeForGamma* fParticleChange;
00127 G4NistManager* fNistManager;
00128 G4IonCoulombCrossSection* ioncross;
00129
00130 const std::vector<G4double>* pCuts;
00131 const G4MaterialCutsCouple* currentCouple;
00132 const G4Material* currentMaterial;
00133 const G4Element* currentElement;
00134 G4int currentMaterialIndex;
00135
00136
00137 G4int heavycorr;
00138
00139 G4double cosThetaMin;
00140 G4double recoilThreshold;
00141
00142
00143
00144 const G4ParticleDefinition* particle;
00145 const G4ParticleDefinition* theProton;
00146 G4double mass;
00147 G4double lowEnergyLimit;
00148
00149
00150 private:
00151
00152 G4bool isInitialised;
00153
00154 };
00155
00156
00157 void G4IonCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
00158 {
00159 if(cup != currentCouple) {
00160 currentCouple = cup;
00161 currentMaterial = cup->GetMaterial();
00162 currentMaterialIndex = currentCouple->GetIndex();
00163
00164 }
00165 }
00166
00167
00168
00169 inline
00170 void G4IonCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
00171 {
00172 if(p != particle) {
00173 particle = p;
00174 mass = particle->GetPDGMass();
00175 ioncross->SetupParticle(p);
00176 }
00177 }
00178
00179
00180 inline void G4IonCoulombScatteringModel::SetRecoilThreshold(G4double eth)
00181 {
00182 recoilThreshold = eth;
00183 }
00184
00185
00186
00187 #endif