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