File indexing completed on 2025-01-18 09:58:44
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 G4NeutronGeneralProcess_h
0048 #define G4NeutronGeneralProcess_h 1
0049
0050 #include "G4HadronicProcess.hh"
0051 #include "globals.hh"
0052 #include "G4HadDataHandler.hh"
0053 #include <vector>
0054
0055 class G4Step;
0056 class G4Track;
0057 class G4ParticleDefinition;
0058 class G4VParticleChange;
0059 class G4VCrossSectionDataSet;
0060 class G4CrossSectionDataStore;
0061
0062
0063
0064 class G4NeutronGeneralProcess : public G4HadronicProcess
0065 {
0066 public:
0067
0068 explicit G4NeutronGeneralProcess(const G4String& pname="NeutronGeneralProc");
0069
0070 ~G4NeutronGeneralProcess() override;
0071
0072 G4bool IsApplicable(const G4ParticleDefinition&) override;
0073
0074 void ProcessDescription(std::ostream& outFile) const override;
0075
0076
0077 void PreparePhysicsTable(const G4ParticleDefinition&) override;
0078
0079
0080 void BuildPhysicsTable(const G4ParticleDefinition&) override;
0081
0082
0083 G4bool StorePhysicsTable(const G4ParticleDefinition* part,
0084 const G4String& directory, G4bool ascii) override;
0085
0086
0087 void StartTracking(G4Track*) override;
0088
0089
0090 G4double PostStepGetPhysicalInteractionLength(
0091 const G4Track& track,
0092 G4double previousStepSize,
0093 G4ForceCondition* condition) override;
0094
0095
0096 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
0097
0098 const G4VProcess* GetCreatorProcess() const override;
0099
0100
0101 const G4String& GetSubProcessName() const;
0102
0103
0104 G4int GetSubProcessSubType() const;
0105
0106 void SetInelasticProcess(G4HadronicProcess*);
0107 void SetElasticProcess(G4HadronicProcess*);
0108 void SetCaptureProcess(G4HadronicProcess*);
0109
0110
0111 G4VCrossSectionDataSet* GetXSection(G4int type);
0112 G4HadronicProcess* GetHadronicProcess(G4int type);
0113
0114 inline const G4VProcess* GetSelectedProcess() const;
0115
0116 inline void SetTimeLimit(G4double val);
0117
0118 inline void SetMinEnergyLimit(G4double val);
0119
0120
0121 G4NeutronGeneralProcess(G4NeutronGeneralProcess &) = delete;
0122 G4NeutronGeneralProcess & operator=
0123 (const G4NeutronGeneralProcess &right) = delete;
0124
0125 protected:
0126
0127 G4double GetMeanFreePath(const G4Track& track, G4double previousStepSize,
0128 G4ForceCondition* condition) override;
0129
0130 inline G4double ComputeGeneralLambda(size_t idxe, size_t idxt);
0131
0132 inline G4double GetProbability(size_t idxt);
0133
0134 inline void SelectedProcess(const G4Step& step, G4HadronicProcess* ptr,
0135 G4CrossSectionDataStore*);
0136
0137 private:
0138
0139
0140 G4double ComputeCrossSection(G4VCrossSectionDataSet*, const G4Material*,
0141 G4double kinEnergy, G4double loge);
0142
0143 G4VCrossSectionDataSet* InitialisationXS(G4HadronicProcess*);
0144
0145
0146 inline void CurrentCrossSection(const G4Track&);
0147
0148 static G4HadDataHandler* theHandler;
0149 static const size_t nTables = 5;
0150 static G4String nameT[nTables];
0151
0152 G4HadronicProcess* fInelasticP = nullptr;
0153 G4HadronicProcess* fElasticP = nullptr;
0154 G4HadronicProcess* fCaptureP = nullptr;
0155 G4HadronicProcess* fSelectedProc = nullptr;
0156
0157 G4VCrossSectionDataSet* fInelasticXS = nullptr;
0158 G4VCrossSectionDataSet* fElasticXS = nullptr;
0159 G4VCrossSectionDataSet* fCaptureXS = nullptr;
0160
0161 G4CrossSectionDataStore* fXSSInelastic = nullptr;
0162 G4CrossSectionDataStore* fXSSElastic = nullptr;
0163 G4CrossSectionDataStore* fXSSCapture = nullptr;
0164 G4CrossSectionDataStore* fCurrentXSS = nullptr;
0165
0166 const G4ParticleDefinition* fNeutron;
0167 const G4Material* fCurrMat = nullptr;
0168
0169 G4double fMinEnergy;
0170 G4double fMiddleEnergy;
0171 G4double fMaxEnergy;
0172 G4double fTimeLimit;
0173 G4double fXSFactorInel = 1.0;
0174 G4double fXSFactorEl = 1.0;
0175 G4double fCurrE = 0.0;
0176 G4double fCurrLogE = 0.0;
0177 G4double fLambda = 0.0;
0178
0179
0180 std::size_t nLowE = 100;
0181 std::size_t nHighE = 10;
0182
0183 std::size_t idxEnergy = 0;
0184 std::size_t matIndex = 0;
0185
0186 G4bool isMaster = true;
0187 std::vector<G4double> fXsec;
0188 };
0189
0190
0191
0192 inline G4double
0193 G4NeutronGeneralProcess::ComputeGeneralLambda(std::size_t idxe, std::size_t idxt)
0194 {
0195 idxEnergy = idxe;
0196 return theHandler->GetVector(idxt, matIndex)
0197 ->LogVectorValue(fCurrE, fCurrLogE);
0198 }
0199
0200
0201
0202 inline G4double G4NeutronGeneralProcess::GetProbability(std::size_t idxt)
0203 {
0204 return theHandler->GetVector(idxt, matIndex)
0205 ->LogVectorValue(fCurrE, fCurrLogE);
0206 }
0207
0208
0209
0210 inline void
0211 G4NeutronGeneralProcess::SelectedProcess(const G4Step& step,
0212 G4HadronicProcess* ptr,
0213 G4CrossSectionDataStore* xs)
0214
0215 {
0216 fSelectedProc = ptr;
0217 fCurrentXSS = xs;
0218 step.GetPostStepPoint()->SetProcessDefinedStep(ptr);
0219 }
0220
0221
0222
0223 inline const G4VProcess* G4NeutronGeneralProcess::GetSelectedProcess() const
0224 {
0225 return fSelectedProc;
0226 }
0227
0228
0229
0230 inline void G4NeutronGeneralProcess::CurrentCrossSection(const G4Track& track)
0231 {
0232 G4double energy = track.GetKineticEnergy();
0233 const G4Material* mat = track.GetMaterial();
0234 if(mat != fCurrMat || energy != fCurrE) {
0235 fCurrMat = mat;
0236 matIndex = mat->GetIndex();
0237 fCurrE = energy;
0238 fCurrLogE = track.GetDynamicParticle()->GetLogKineticEnergy();
0239 fLambda = (energy <= fMiddleEnergy) ? ComputeGeneralLambda(0, 0)
0240 : ComputeGeneralLambda(1, 3);
0241 currentInteractionLength = 1.0/fLambda;
0242 }
0243 }
0244
0245
0246
0247 inline void G4NeutronGeneralProcess::SetTimeLimit(G4double val)
0248 {
0249 fTimeLimit = val;
0250 }
0251
0252
0253
0254 inline void G4NeutronGeneralProcess::SetMinEnergyLimit(G4double val)
0255 {
0256 fMinEnergy = val;
0257 }
0258
0259
0260
0261 #endif