File indexing completed on 2025-01-18 09:58:31
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057 #ifndef G4IonCoulombScatteringModel_h
0058 #define G4IonCoulombScatteringModel_h 1
0059
0060 #include "G4VEmModel.hh"
0061 #include "globals.hh"
0062 #include "G4NistManager.hh"
0063 #include "G4IonCoulombCrossSection.hh"
0064
0065 #include <vector>
0066
0067 class G4ParticleChangeForGamma;
0068 class G4ParticleDefinition;
0069 class G4IonTable;
0070
0071 class G4IonCoulombScatteringModel : public G4VEmModel
0072 {
0073 public:
0074
0075 explicit G4IonCoulombScatteringModel(const G4String& nam =
0076 "IonCoulombScattering");
0077
0078 ~G4IonCoulombScatteringModel() override;
0079
0080 void Initialise(const G4ParticleDefinition*, const G4DataVector&) final;
0081
0082 G4double ComputeCrossSectionPerAtom(
0083 const G4ParticleDefinition*,
0084 G4double kinEnergy,
0085 G4double Z,
0086 G4double A,
0087 G4double cut,
0088 G4double emax) final;
0089
0090 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0091 const G4MaterialCutsCouple*,
0092 const G4DynamicParticle*,
0093 G4double tmin,
0094 G4double maxEnergy) final;
0095
0096
0097
0098 inline void SetRecoilThreshold(G4double eth);
0099
0100 inline void SetHeavyIonCorr(G4int b);
0101
0102 inline G4int GetHeavyIonCorr();
0103
0104
0105 G4IonCoulombScatteringModel & operator=
0106 (const G4IonCoulombScatteringModel &right) = delete;
0107 G4IonCoulombScatteringModel(const G4IonCoulombScatteringModel&) = delete;
0108
0109 private:
0110
0111 inline void DefineMaterial(const G4MaterialCutsCouple*);
0112
0113 inline void SetupParticle(const G4ParticleDefinition*);
0114
0115 G4IonTable* theIonTable;
0116 G4ParticleChangeForGamma* fParticleChange;
0117 G4NistManager* fNistManager;
0118 G4IonCoulombCrossSection* ioncross;
0119
0120 const std::vector<G4double>* pCuts;
0121 const G4MaterialCutsCouple* currentCouple;
0122 const G4Material* currentMaterial;
0123 const G4Element* currentElement;
0124 G4int currentMaterialIndex;
0125
0126 G4int heavycorr;
0127
0128 G4double cosThetaMin;
0129 G4double recoilThreshold;
0130
0131
0132 const G4ParticleDefinition* particle;
0133 const G4ParticleDefinition* theProton;
0134 G4double mass;
0135
0136 };
0137
0138
0139
0140 inline void
0141 G4IonCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
0142 {
0143 if(cup != currentCouple) {
0144 currentCouple = cup;
0145 currentMaterial = cup->GetMaterial();
0146 currentMaterialIndex = currentCouple->GetIndex();
0147 }
0148 }
0149
0150
0151
0152 inline
0153 void G4IonCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
0154 {
0155 if(p != particle) {
0156 particle = p;
0157 mass = particle->GetPDGMass();
0158 ioncross->SetupParticle(p);
0159 }
0160 }
0161
0162
0163
0164 inline void G4IonCoulombScatteringModel::SetRecoilThreshold(G4double eth)
0165 {
0166 recoilThreshold = eth;
0167 }
0168
0169
0170
0171 inline void G4IonCoulombScatteringModel::SetHeavyIonCorr(G4int b)
0172 {
0173 heavycorr = b;
0174 }
0175
0176
0177
0178 inline G4int G4IonCoulombScatteringModel::GetHeavyIonCorr()
0179 {
0180 return heavycorr;
0181 }
0182
0183
0184
0185 #endif