File indexing completed on 2025-01-18 09:58:48
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 #ifndef G4PAIModel_h
0054 #define G4PAIModel_h 1
0055
0056 #include <CLHEP/Units/PhysicalConstants.h>
0057
0058 #include "G4VEmModel.hh"
0059 #include "G4VEmFluctuationModel.hh"
0060 #include "globals.hh"
0061 #include <vector>
0062
0063 class G4Region;
0064 class G4MaterialCutsCouple;
0065 class G4ParticleChangeForLoss;
0066 class G4PAIModelData;
0067
0068 class G4PAIModel final : public G4VEmModel, public G4VEmFluctuationModel
0069 {
0070
0071 public:
0072
0073 explicit G4PAIModel(const G4ParticleDefinition* p = nullptr,
0074 const G4String& nam = "PAI");
0075
0076 ~G4PAIModel() final;
0077
0078 void Initialise(const G4ParticleDefinition*, const G4DataVector&) final;
0079
0080 void InitialiseLocal(const G4ParticleDefinition*,
0081 G4VEmModel* masterModel) final;
0082
0083 G4double MinEnergyCut(const G4ParticleDefinition*,
0084 const G4MaterialCutsCouple* couple) final;
0085
0086 G4double ComputeDEDXPerVolume(const G4Material*,
0087 const G4ParticleDefinition*,
0088 G4double kineticEnergy,
0089 G4double cutEnergy) final;
0090
0091 G4double CrossSectionPerVolume(const G4Material*,
0092 const G4ParticleDefinition*,
0093 G4double kineticEnergy,
0094 G4double cutEnergy,
0095 G4double maxEnergy) final;
0096
0097 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0098 const G4MaterialCutsCouple*,
0099 const G4DynamicParticle*,
0100 G4double tmin,
0101 G4double maxEnergy) final;
0102
0103 G4double SampleFluctuations(const G4MaterialCutsCouple*,
0104 const G4DynamicParticle*,
0105 const G4double, const G4double,
0106 const G4double, const G4double) final;
0107
0108 G4double Dispersion(const G4Material*, const G4DynamicParticle*,
0109 const G4double, const G4double,
0110 const G4double) final;
0111
0112 void DefineForRegion(const G4Region* r) final;
0113
0114 inline G4PAIModelData* GetPAIModelData();
0115
0116 inline const std::vector<const G4MaterialCutsCouple*>& GetVectorOfCouples();
0117
0118 inline G4double ComputeMaxEnergy(G4double scaledEnergy);
0119
0120 inline void SetVerboseLevel(G4int verbose);
0121
0122
0123 G4PAIModel & operator=(const G4PAIModel &right) = delete;
0124 G4PAIModel(const G4PAIModel&) = delete;
0125
0126 protected:
0127
0128 G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
0129 G4double kinEnergy) final;
0130
0131 private:
0132
0133 inline G4int FindCoupleIndex(const G4MaterialCutsCouple*);
0134
0135 inline void SetParticle(const G4ParticleDefinition* p);
0136
0137 G4int fVerbose;
0138
0139 G4PAIModelData* fModelData;
0140
0141 std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
0142 std::vector<const G4Region*> fPAIRegionVector;
0143
0144 const G4ParticleDefinition* fParticle;
0145 const G4ParticleDefinition* fElectron;
0146 const G4ParticleDefinition* fPositron;
0147 G4ParticleChangeForLoss* fParticleChange;
0148
0149 G4double fMass;
0150 G4double fRatio;
0151 G4double fChargeSquare;
0152 G4double fLowestTcut;
0153 };
0154
0155 inline G4PAIModelData* G4PAIModel::GetPAIModelData()
0156 {
0157 return fModelData;
0158 }
0159
0160 inline const std::vector<const G4MaterialCutsCouple*>&
0161 G4PAIModel::GetVectorOfCouples()
0162 {
0163 return fMaterialCutsCoupleVector;
0164 }
0165
0166 inline G4double G4PAIModel::ComputeMaxEnergy(G4double scaledEnergy)
0167 {
0168 return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
0169 }
0170
0171 inline void G4PAIModel::SetVerboseLevel(G4int verbose)
0172 {
0173 fVerbose=verbose;
0174 }
0175
0176 inline G4int G4PAIModel::FindCoupleIndex(const G4MaterialCutsCouple* couple)
0177 {
0178 G4int idx = -1;
0179 G4int jMatMax = (G4int)fMaterialCutsCoupleVector.size();
0180 for(G4int jMat = 0;jMat < jMatMax; ++jMat) {
0181 if(couple == fMaterialCutsCoupleVector[jMat]) {
0182 idx = jMat;
0183 break;
0184 }
0185 }
0186 return idx;
0187 }
0188
0189 inline void G4PAIModel::SetParticle(const G4ParticleDefinition* p)
0190 {
0191 if(fParticle != p) {
0192 fParticle = p;
0193 fMass = fParticle->GetPDGMass();
0194 fRatio = CLHEP::proton_mass_c2/fMass;
0195 G4double q = fParticle->GetPDGCharge()/CLHEP::eplus;
0196 fChargeSquare = q*q;
0197 }
0198 }
0199
0200 #endif