Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:44

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // -------------------------------------------------------------------
0028 //
0029 // GEANT4 Class header file
0030 //
0031 //
0032 // File name:     G4NeutronGeneralProcess
0033 //
0034 // Author:        Vladimir Ivanchenko
0035 //
0036 // Creation date: 08.08.2022
0037 //
0038 // Modifications:
0039 //
0040 // Class Description:
0041 //
0042 // It is the neutron super process
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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   // Initialise for build of tables
0077   void PreparePhysicsTable(const G4ParticleDefinition&) override;
0078 
0079   // Build physics table during initialisation
0080   void BuildPhysicsTable(const G4ParticleDefinition&) override;
0081 
0082   // Store internal tables after initialisation
0083   G4bool StorePhysicsTable(const G4ParticleDefinition* part,
0084                            const G4String& directory, G4bool ascii) override;
0085 
0086   // Called before tracking of each new G4Track
0087   void StartTracking(G4Track*) override;
0088 
0089   // implementation of virtual method, specific for G4NeutronGeneralProcess
0090   G4double PostStepGetPhysicalInteractionLength(
0091                              const G4Track& track,
0092                              G4double previousStepSize,
0093                              G4ForceCondition* condition) override;
0094 
0095   // implementation of virtual method, specific for G4NeutronGeneralProcess
0096   G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
0097 
0098   const G4VProcess* GetCreatorProcess() const override;
0099 
0100   // Temporary method
0101   const G4String& GetSubProcessName() const;
0102 
0103   // Temporary method
0104   G4int GetSubProcessSubType() const;
0105 
0106   void SetInelasticProcess(G4HadronicProcess*);
0107   void SetElasticProcess(G4HadronicProcess*);
0108   void SetCaptureProcess(G4HadronicProcess*);
0109 
0110   // access methods to cross sections and processes
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   // hide copy constructor and assignment operator
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   // partial cross section
0140   G4double ComputeCrossSection(G4VCrossSectionDataSet*, const G4Material*,
0141                                G4double kinEnergy, G4double loge);
0142 
0143   G4VCrossSectionDataSet* InitialisationXS(G4HadronicProcess*);
0144 
0145   // total cross section
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   // number of bins per decade
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0201 
0202 inline G4double G4NeutronGeneralProcess::GetProbability(std::size_t idxt)
0203 {
0204   return theHandler->GetVector(idxt, matIndex)
0205     ->LogVectorValue(fCurrE, fCurrLogE);
0206 }
0207 
0208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0222 
0223 inline const G4VProcess* G4NeutronGeneralProcess::GetSelectedProcess() const
0224 {
0225   return fSelectedProc;
0226 }
0227 
0228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0246 
0247 inline void G4NeutronGeneralProcess::SetTimeLimit(G4double val)
0248 {
0249   fTimeLimit = val;
0250 }
0251 
0252 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0253 
0254 inline void G4NeutronGeneralProcess::SetMinEnergyLimit(G4double val)
0255 {
0256   fMinEnergy = val;
0257 }
0258 
0259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0260 
0261 #endif