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
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043 #ifndef G4DNACPA100IonisationModel_h
0044 #define G4DNACPA100IonisationModel_h 1
0045
0046 #include "G4DNACPA100IonisationStructure.hh"
0047 #include "G4DNACrossSectionDataSet.hh"
0048 #include "G4Electron.hh"
0049 #include "G4LogLogInterpolation.hh"
0050 #include "G4NistManager.hh"
0051 #include "G4ParticleChangeForGamma.hh"
0052 #include "G4ProductionCutsTable.hh"
0053 #include "G4VAtomDeexcitation.hh"
0054 #include "G4VDNAModel.hh"
0055
0056 class G4DNACPA100IonisationModel : public G4VDNAModel
0057 {
0058 using TriDimensionMap =
0059 std::map<std::size_t,
0060 std::map<const G4ParticleDefinition*,
0061 std::map<G4int, std::map<G4double, std::map<G4double, G4double>>>>>;
0062 using VecMap =
0063 std::map<std::size_t,
0064 std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>;
0065 using VecMapWithShell =
0066 std::map<std::size_t,
0067 std::map<const G4ParticleDefinition*,
0068 std::map<G4double, std::map<G4double, std::vector<G4double>>>>>;
0069 using PartKineticInMat =
0070 const std::tuple<std::size_t , G4double , G4int >&;
0071
0072 public:
0073 explicit G4DNACPA100IonisationModel(const G4ParticleDefinition* p = nullptr,
0074 const G4String& nam = "DNACPA100IonisationModel");
0075
0076 ~G4DNACPA100IonisationModel() override = default;
0077
0078 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0079
0080 G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p,
0081 G4double ekin, G4double emin, G4double emax) override;
0082
0083 void SampleSecondaries(std::vector<G4DynamicParticle*>*, const G4MaterialCutsCouple*,
0084 const G4DynamicParticle*, G4double tmin, G4double maxEnergy) override;
0085
0086 G4double DifferentialCrossSection(PartKineticInMat info, const G4double& energyTransfer);
0087
0088
0089
0090
0091
0092 inline void SelectFasterComputation(G4bool input);
0093
0094 inline void SelectUseDcs(G4bool input);
0095
0096 inline void SelectStationary(G4bool input);
0097
0098 G4DNACPA100IonisationModel& operator=(const G4DNACPA100IonisationModel& right) = delete;
0099 G4DNACPA100IonisationModel(const G4DNACPA100IonisationModel&) = delete;
0100 void ReadDiffCSFile(const std::size_t& materialID, const G4ParticleDefinition* p,
0101 const G4String& file, const G4double& scaleFactor) override;
0102
0103 protected:
0104 G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr;
0105
0106 private:
0107 G4bool statCode = false;
0108 G4bool fasterCode = true;
0109 G4bool useDcs = false;
0110
0111
0112
0113
0114 G4VAtomDeexcitation* fAtomDeexcitation = nullptr;
0115
0116 G4bool isInitialised = false;
0117 G4int verboseLevel = 0;
0118
0119 G4DNACPA100IonisationStructure iStructure;
0120
0121 G4double RandomizeEjectedElectronEnergy(PartKineticInMat info);
0122
0123 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(PartKineticInMat info);
0124
0125 G4double RandomizeEjectedElectronEnergyFromanalytical(PartKineticInMat info);
0126
0127 G4double RandomTransferedEnergy(PartKineticInMat info);
0128
0129 void RandomizeEjectedElectronDirection(G4ParticleDefinition* aParticleDefinition,
0130 G4double incomingParticleEnergy,
0131 G4double outgoingParticleEnergy, G4double& cosTheta,
0132 G4double& phi);
0133
0134 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
0135
0136 G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11,
0137 G4double x12, G4double x21, G4double x22, G4double t1, G4double t2,
0138 G4double t, G4double e);
0139
0140 TriDimensionMap diffCrossSectionData, fEnergySecondaryData;
0141 std::map<std::size_t, std::map<const G4ParticleDefinition*, std::vector<G4double>>>
0142 fTMapWithVec;
0143 VecMap fEMapWithVector;
0144 VecMapWithShell fProbaShellMap;
0145
0146 const G4Material* fpGuanine = nullptr;
0147 const G4Material* fpG4_WATER = nullptr;
0148 const G4Material* fpDeoxyribose = nullptr;
0149 const G4Material* fpCytosine = nullptr;
0150 const G4Material* fpThymine = nullptr;
0151 const G4Material* fpAdenine = nullptr;
0152 const G4Material* fpPhosphate = nullptr;
0153 const G4ParticleDefinition* fpParticle = nullptr;
0154 G4DNACPA100IonisationModel* fpModelData = nullptr;
0155 };
0156
0157
0158
0159 inline void G4DNACPA100IonisationModel::SelectFasterComputation(G4bool input)
0160 {
0161 fasterCode = input;
0162 }
0163
0164
0165
0166 inline void G4DNACPA100IonisationModel::SelectUseDcs(G4bool input)
0167 {
0168 useDcs = input;
0169 }
0170
0171
0172
0173
0174
0175 inline void G4DNACPA100IonisationModel::SelectStationary(G4bool input)
0176 {
0177 statCode = input;
0178 }
0179
0180
0181
0182 #endif