File indexing completed on 2025-01-18 09:58:19
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 #ifndef G4GammaGeneralProcess_h
0048 #define G4GammaGeneralProcess_h 1
0049
0050 #include "G4VEmProcess.hh"
0051 #include "globals.hh"
0052 #include "G4EmDataHandler.hh"
0053
0054 class G4Step;
0055 class G4Track;
0056 class G4ParticleDefinition;
0057 class G4VParticleChange;
0058 class G4GammaConversionToMuons;
0059 class G4HadronicProcess;
0060 class G4MaterialCutsCouple;
0061
0062
0063
0064 class G4GammaGeneralProcess : public G4VEmProcess
0065 {
0066 public:
0067
0068 explicit G4GammaGeneralProcess(const G4String& pname="GammaGeneralProc");
0069
0070 ~G4GammaGeneralProcess() override;
0071
0072 G4bool IsApplicable(const G4ParticleDefinition&) override;
0073
0074 void AddEmProcess(G4VEmProcess*);
0075
0076 void AddMMProcess(G4GammaConversionToMuons*);
0077
0078 void AddHadProcess(G4HadronicProcess*);
0079
0080 void ProcessDescription(std::ostream& outFile) const override;
0081
0082 protected:
0083
0084 void InitialiseProcess(const G4ParticleDefinition*) override;
0085
0086 public:
0087
0088
0089 void PreparePhysicsTable(const G4ParticleDefinition&) override;
0090
0091
0092 void BuildPhysicsTable(const G4ParticleDefinition&) override;
0093
0094
0095 void StartTracking(G4Track*) override;
0096
0097
0098 G4double PostStepGetPhysicalInteractionLength(
0099 const G4Track& track,
0100 G4double previousStepSize,
0101 G4ForceCondition* condition) override;
0102
0103
0104 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
0105
0106
0107
0108 G4bool StorePhysicsTable(const G4ParticleDefinition*,
0109 const G4String& directory,
0110 G4bool ascii = false) override;
0111
0112
0113
0114
0115
0116
0117 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
0118 const G4String& directory,
0119 G4bool ascii) override;
0120
0121
0122 const G4VProcess* GetCreatorProcess() const override;
0123 inline const G4VProcess* GetSelectedProcess() const;
0124
0125
0126 const G4String& GetSubProcessName() const;
0127 G4int GetSubProcessSubType() const;
0128
0129 G4VEmProcess* GetEmProcess(const G4String& name) override;
0130
0131 inline G4HadronicProcess* GetGammaNuclear() const;
0132
0133
0134 G4GammaGeneralProcess(G4GammaGeneralProcess &) = delete;
0135 G4GammaGeneralProcess & operator=
0136 (const G4GammaGeneralProcess &right) = delete;
0137
0138 protected:
0139
0140 G4double GetMeanFreePath(const G4Track& track, G4double previousStepSize,
0141 G4ForceCondition* condition) override;
0142
0143 inline G4double ComputeGeneralLambda(size_t idxe, size_t idxt);
0144
0145 inline G4double GetProbability(size_t idxt);
0146
0147 inline void SelectedProcess(const G4Step& step, G4VProcess* ptr);
0148
0149 inline void SelectEmProcess(const G4Step&, G4VEmProcess*);
0150
0151 void SelectHadProcess(const G4Track&, const G4Step&, G4HadronicProcess*);
0152
0153
0154 G4double TotalCrossSectionPerVolume();
0155
0156 private:
0157
0158 G4bool RetrieveTable(G4VEmProcess*, const G4String& directory,
0159 G4bool ascii);
0160
0161 protected:
0162
0163 G4HadronicProcess* theGammaNuclear = nullptr;
0164 G4VProcess* selectedProc = nullptr;
0165
0166 G4double preStepLogE = 1.0;
0167 G4double factor = 1.0;
0168
0169
0170 private:
0171 static G4EmDataHandler* theHandler;
0172 static const size_t nTables = 15;
0173 static G4bool theT[nTables];
0174 static G4String nameT[nTables];
0175
0176 G4VEmProcess* thePhotoElectric = nullptr;
0177 G4VEmProcess* theCompton = nullptr;
0178 G4VEmProcess* theConversionEE = nullptr;
0179 G4VEmProcess* theRayleigh = nullptr;
0180 G4GammaConversionToMuons* theConversionMM = nullptr;
0181
0182 G4double minPEEnergy;
0183 G4double minEEEnergy;
0184 G4double minMMEnergy;
0185 G4double peLambda = 0.0;
0186
0187 size_t nLowE = 40;
0188 size_t nHighE = 50;
0189 size_t idxEnergy = 0;
0190 };
0191
0192
0193
0194 inline G4double
0195 G4GammaGeneralProcess::ComputeGeneralLambda(size_t idxe, size_t idxt)
0196 {
0197 idxEnergy = idxe;
0198 return factor*theHandler->GetVector(idxt, basedCoupleIndex)
0199 ->LogVectorValue(preStepKinEnergy, preStepLogE);
0200 }
0201
0202
0203
0204 inline G4double G4GammaGeneralProcess::GetProbability(size_t idxt)
0205 {
0206 return theHandler->GetVector(idxt, basedCoupleIndex)
0207 ->LogVectorValue(preStepKinEnergy, preStepLogE);
0208 }
0209
0210
0211
0212 inline void
0213 G4GammaGeneralProcess::SelectedProcess(const G4Step& step, G4VProcess* ptr)
0214 {
0215 selectedProc = ptr;
0216 step.GetPostStepPoint()->SetProcessDefinedStep(ptr);
0217 }
0218
0219
0220
0221 inline void G4GammaGeneralProcess::SelectEmProcess(const G4Step& step, G4VEmProcess* proc)
0222 {
0223 proc->CurrentSetup(currentCouple, preStepKinEnergy);
0224 SelectedProcess(step, proc);
0225 }
0226
0227
0228
0229 inline const G4VProcess* G4GammaGeneralProcess::GetSelectedProcess() const
0230 {
0231 return selectedProc;
0232 }
0233
0234
0235
0236
0237 inline G4HadronicProcess* G4GammaGeneralProcess::GetGammaNuclear() const
0238 {
0239 return theGammaNuclear;
0240 }
0241
0242
0243
0244 #endif