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 G4DNABornIonisationModel1_h
0029 #define G4DNABornIonisationModel1_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 G4DNABornIonisationModel1 : public G4VEmModel
0048 {
0049
0050 public:
0051
0052 G4DNABornIonisationModel1(const G4ParticleDefinition* p = nullptr,
0053 const G4String& nam = "DNABornIonisationModel");
0054
0055 ~G4DNABornIonisationModel1() override;
0056
0057 G4DNABornIonisationModel1 & operator=(const G4DNABornIonisationModel1 &right) = delete;
0058 G4DNABornIonisationModel1(const G4DNABornIonisationModel1&) = 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,
0082 G4double incomingParticleEnergy, G4int shell, G4double random) ;
0083
0084 inline void SelectFasterComputation(G4bool input);
0085
0086 inline void SelectStationary(G4bool input);
0087
0088 inline void SelectSPScaling(G4bool input);
0089
0090 protected:
0091
0092 G4ParticleChangeForGamma* fParticleChangeForGamma;
0093
0094 private:
0095
0096 G4bool fasterCode;
0097 G4bool statCode;
0098 G4bool spScaling;
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
0110
0111
0112
0113 G4bool isInitialised{false};
0114 G4int verboseLevel;
0115
0116
0117
0118 using MapFile = std::map<G4String, G4String, std::less<G4String>>;
0119 MapFile tableFile;
0120
0121 using MapData = std::map<G4String, G4DNACrossSectionDataSet *, std::less<G4String>>;
0122 MapData tableData;
0123
0124
0125
0126 G4DNAWaterIonisationStructure waterStructure;
0127
0128 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
0129
0130 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
0131
0132 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
0133
0134 G4double QuadInterpolator( G4double e11,
0135 G4double e12,
0136 G4double e21,
0137 G4double e22,
0138 G4double x11,
0139 G4double x12,
0140 G4double x21,
0141 G4double x22,
0142 G4double t1,
0143 G4double t2,
0144 G4double t,
0145 G4double e);
0146
0147 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
0148
0149 TriDimensionMap eDiffCrossSectionData[6];
0150 TriDimensionMap eNrjTransfData[6];
0151
0152 TriDimensionMap pDiffCrossSectionData[6];
0153 TriDimensionMap pNrjTransfData[6];
0154
0155 std::vector<G4double> eTdummyVec;
0156 std::vector<G4double> pTdummyVec;
0157
0158 using VecMap = std::map<G4double, std::vector<G4double>>;
0159
0160 VecMap eVecm;
0161 VecMap pVecm;
0162
0163 VecMap eProbaShellMap[6];
0164 VecMap pProbaShellMap[6];
0165
0166
0167
0168 G4int RandomSelect(G4double energy,const G4String& particle );
0169
0170 };
0171
0172 inline void G4DNABornIonisationModel1::SelectFasterComputation (G4bool input)
0173 {
0174 fasterCode = input;
0175 }
0176
0177
0178
0179 inline void G4DNABornIonisationModel1::SelectStationary (G4bool input)
0180 {
0181 statCode = input;
0182 }
0183
0184
0185
0186 inline void G4DNABornIonisationModel1::SelectSPScaling (G4bool input)
0187 {
0188 spScaling = input;
0189 }
0190
0191
0192
0193 #endif