File indexing completed on 2025-01-18 09:58:07
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 #ifndef G4DNABornIonisationModel2_h
0029 #define G4DNABornIonisationModel2_h 1
0030
0031 #include "G4VEmModel.hh"
0032 #include "G4ParticleChangeForGamma.hh"
0033 #include "G4ProductionCutsTable.hh"
0034
0035 #include "G4DNACrossSectionDataSet.hh"
0036 #include "G4Electron.hh"
0037 #include "G4Proton.hh"
0038 #include "G4DNAGenericIonsManager.hh"
0039
0040 #include "G4LogLogInterpolation.hh"
0041
0042 #include "G4DNAWaterIonisationStructure.hh"
0043 #include "G4VAtomDeexcitation.hh"
0044 #include "G4NistManager.hh"
0045
0046
0047 class G4DNABornIonisationModel2 : public G4VEmModel
0048 {
0049
0050 public:
0051
0052 G4DNABornIonisationModel2(const G4ParticleDefinition* p = nullptr,
0053 const G4String& nam = "DNABornIonisationModel");
0054
0055 ~G4DNABornIonisationModel2() override;
0056
0057 G4DNABornIonisationModel2 & operator=(const G4DNABornIonisationModel2 &right) = delete;
0058 G4DNABornIonisationModel2(const G4DNABornIonisationModel2&) = delete;
0059
0060 void Initialise(const G4ParticleDefinition*, const G4DataVector& = *(new G4DataVector())) override;
0061
0062 G4double CrossSectionPerVolume( const G4Material* material,
0063 const G4ParticleDefinition* p,
0064 G4double ekin,
0065 G4double emin,
0066 G4double emax) override;
0067
0068 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0069 const G4MaterialCutsCouple*,
0070 const G4DynamicParticle*,
0071 G4double tmin,
0072 G4double maxEnergy) override;
0073
0074 G4double GetPartialCrossSection(const G4Material*,
0075 G4int ,
0076 const G4ParticleDefinition*,
0077 G4double ) override;
0078
0079 G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
0080
0081 G4double TransferedEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random) ;
0082
0083 inline void SelectFasterComputation(G4bool input);
0084
0085 inline void SelectStationary(G4bool input);
0086
0087 inline void SelectSPScaling(G4bool input);
0088
0089 protected:
0090
0091 G4ParticleChangeForGamma* fParticleChangeForGamma;
0092
0093 private:
0094
0095 G4bool fasterCode;
0096 G4bool statCode;
0097 G4bool spScaling;
0098
0099
0100 const std::vector<G4double>* fpMolWaterDensity;
0101
0102
0103 G4VAtomDeexcitation* fAtomDeexcitation;
0104
0105 G4double fLowEnergyLimit;
0106 G4double fHighEnergyLimit;
0107
0108 const G4ParticleDefinition* fParticleDef;
0109
0110 G4bool isInitialised{false};
0111 G4int verboseLevel;
0112
0113
0114 G4String fTableFile;
0115 G4DNACrossSectionDataSet* fTableData;
0116
0117
0118
0119 G4DNAWaterIonisationStructure waterStructure;
0120
0121 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
0122
0123 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
0124
0125 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
0126
0127 G4double QuadInterpolator( G4double e11,
0128 G4double e12,
0129 G4double e21,
0130 G4double e22,
0131 G4double x11,
0132 G4double x12,
0133 G4double x21,
0134 G4double x22,
0135 G4double t1,
0136 G4double t2,
0137 G4double t,
0138 G4double e);
0139
0140 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
0141
0142 TriDimensionMap fDiffCrossSectionData[6];
0143 TriDimensionMap fNrjTransfData[6];
0144
0145 std::vector<G4double> fTdummyVec;
0146
0147 using VecMap = std::map<G4double, std::vector<G4double>>;
0148
0149 VecMap fVecm;
0150 VecMap fProbaShellMap[6];
0151
0152
0153
0154 G4int RandomSelect(G4double energy);
0155
0156 };
0157
0158
0159
0160 inline void G4DNABornIonisationModel2::SelectFasterComputation (G4bool input)
0161 {
0162 fasterCode = input;
0163 }
0164
0165
0166
0167 inline void G4DNABornIonisationModel2::SelectStationary (G4bool input)
0168 {
0169 statCode = input;
0170 }
0171
0172
0173
0174 inline void G4DNABornIonisationModel2::SelectSPScaling (G4bool input)
0175 {
0176 spScaling = input;
0177 }
0178
0179
0180
0181 #endif