Back to home page

EIC code displayed by LXR

 
 

    


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

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 // G4HadronicProcess
0032 //
0033 // This is the top level Hadronic Process class
0034 // The inelastic, elastic, capture, and fission processes
0035 // should derive from this class
0036 //
0037 // original by H.P.Wellisch
0038 // J.L. Chuma, TRIUMF, 10-Mar-1997
0039 // Last modified: 04-Apr-1997
0040 // 19-May-2008 V.Ivanchenko cleanup and added comments
0041 // 05-Jul-2010 V.Ivanchenko cleanup commented lines
0042 // 28-Jul-2012 M.Maire add function GetTargetDefinition() 
0043 // 14-Sep-2012 Inherit from RestDiscrete, use subtype code (now in ctor) to
0044 //      configure base-class
0045 // 28-Sep-2012 M. Kelsey -- Undo inheritance change, keep new ctor
0046 
0047 #ifndef G4HadronicProcess_h
0048 #define G4HadronicProcess_h 1
0049  
0050 #include "globals.hh"
0051 #include "G4VDiscreteProcess.hh"
0052 #include "G4EnergyRangeManager.hh"
0053 #include "G4Nucleus.hh" 
0054 #include "G4ReactionProduct.hh"
0055 #include "G4HadronicProcessType.hh"
0056 #include "G4CrossSectionDataStore.hh"
0057 #include "G4Material.hh"
0058 #include "G4DynamicParticle.hh"
0059 #include "G4ThreeVector.hh"
0060 #include "G4HadXSTypes.hh"
0061 #include <vector>
0062 
0063 class G4Track;
0064 class G4Step;
0065 class G4Element;
0066 class G4ParticleChange;
0067 class G4HadronicInteraction;
0068 class G4HadronicProcessStore;
0069 class G4VCrossSectionDataSet;
0070 class G4VLeadingParticleBiasing;
0071 class G4ParticleDefinition;
0072 
0073 class G4HadronicProcess : public G4VDiscreteProcess
0074 {
0075 public:
0076   G4HadronicProcess(const G4String& processName="Hadronic",
0077             G4ProcessType procType=fHadronic);    
0078 
0079   // Preferred signature for subclasses, specifying their subtype here
0080   G4HadronicProcess(const G4String& processName, 
0081             G4HadronicProcessType subType);    
0082 
0083   ~G4HadronicProcess() override;
0084 
0085   // register generator of secondaries
0086   void RegisterMe(G4HadronicInteraction* a);
0087 
0088   // get cross section per element
0089   G4double GetElementCrossSection(const G4DynamicParticle * part, 
0090                   const G4Element * elm, 
0091                   const G4Material* mat = nullptr);
0092 
0093   // obsolete method to get cross section per element
0094   inline
0095   G4double GetMicroscopicCrossSection(const G4DynamicParticle * part, 
0096                       const G4Element * elm, 
0097                       const G4Material* mat = nullptr);
0098 
0099   // initialisation for a new track
0100   void StartTracking(G4Track* track) override;
0101 
0102   // compute step limit
0103   G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
0104            G4double, G4ForceCondition*) override;
0105 
0106   // generic PostStepDoIt recommended for all derived classes
0107   G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
0108                   const G4Step& aStep) override;
0109 
0110   // initialisation of physics tables and G4HadronicProcessStore
0111   void PreparePhysicsTable(const G4ParticleDefinition&) override;
0112 
0113   // build physics tables and print out the configuration of the process
0114   void BuildPhysicsTable(const G4ParticleDefinition&) override;
0115 
0116   // dump physics tables 
0117   void DumpPhysicsTable(const G4ParticleDefinition& p);
0118 
0119   // add cross section data set
0120   void AddDataSet(G4VCrossSectionDataSet * aDataSet);
0121 
0122   // access to the list of hadronic interactions
0123   std::vector<G4HadronicInteraction*>& GetHadronicInteractionList();
0124 
0125   // access to an hadronic interaction by name
0126   G4HadronicInteraction* GetHadronicModel(const G4String&);
0127 
0128   // access to the chosen generator
0129   inline G4HadronicInteraction* GetHadronicInteraction() const;
0130   
0131   // get inverse cross section per volume
0132   G4double GetMeanFreePath(const G4Track &aTrack, G4double, 
0133                G4ForceCondition *) override;
0134 
0135   // access to the target nucleus
0136   inline const G4Nucleus* GetTargetNucleus() const;
0137 
0138   inline G4Nucleus* GetTargetNucleusPointer();
0139   
0140   inline const G4Isotope* GetTargetIsotope();
0141 
0142   // methods needed for implementation of integral XS
0143   G4double ComputeCrossSection(const G4ParticleDefinition*,
0144                                const G4Material*,
0145                                const G4double kinEnergy);
0146 
0147   inline G4HadXSType CrossSectionType() const;
0148   inline void SetCrossSectionType(G4HadXSType val);
0149 
0150   void ProcessDescription(std::ostream& outFile) const override;
0151  
0152   // scale cross section
0153   void BiasCrossSectionByFactor(G4double aScale);
0154   void MultiplyCrossSectionBy(G4double factor);
0155   inline G4double CrossSectionFactor() const;
0156 
0157   // Integral option 
0158   inline void SetIntegral(G4bool val);
0159 
0160   // Energy-momentum non-conservation limits and reporting
0161   inline void SetEpReportLevel(G4int level);
0162   inline void SetEnergyMomentumCheckLevels(G4double relativeLevel,
0163                                            G4double absoluteLevel);
0164   inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const;
0165 
0166   // access to the cross section data store
0167   inline G4CrossSectionDataStore* GetCrossSectionDataStore();
0168 
0169   // access to the data for integral XS method
0170   inline std::vector<G4TwoPeaksHadXS*>* TwoPeaksXS() const;
0171   inline std::vector<G4double>* EnergyOfCrossSectionMax() const;
0172 
0173   // hide assignment operator as private 
0174   G4HadronicProcess& operator=(const G4HadronicProcess& right) = delete;
0175   G4HadronicProcess(const G4HadronicProcess&) = delete;
0176 
0177 protected:
0178 
0179   // generic method to choose secondary generator 
0180   // recommended for all derived classes
0181   inline G4HadronicInteraction* ChooseHadronicInteraction(
0182       const G4HadProjectile & aHadProjectile, G4Nucleus& aTargetNucleus,
0183       const G4Material* aMaterial, const G4Element* anElement);
0184               
0185   // access to the cross section data set
0186   inline G4double GetLastCrossSection();
0187 
0188   // fill result
0189   void FillResult(G4HadFinalState* aR, const G4Track& aT);
0190 
0191   void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
0192 
0193   // Check the result for catastrophic energy non-conservation
0194   G4HadFinalState* CheckResult(const G4HadProjectile& thePro,
0195                    const G4Nucleus& targetNucleus, 
0196                    G4HadFinalState* result);
0197 
0198   // Check 4-momentum balance
0199   void CheckEnergyMomentumConservation(const G4Track&, const G4Nucleus&);
0200 
0201 private:
0202 
0203   void InitialiseLocal();
0204   void UpdateCrossSectionAndMFP(const G4double kinEnergy);
0205   void RecomputeXSandMFP(const G4double kinEnergy);
0206 
0207   inline void DefineXSandMFP();
0208   inline void ComputeXSandMFP();
0209 
0210   G4double XBiasSurvivalProbability();
0211   G4double XBiasSecondaryWeight();
0212 
0213 protected:
0214 
0215   G4HadProjectile thePro;
0216 
0217   G4ParticleChange* theTotalResult;
0218   G4CrossSectionDataStore* theCrossSectionDataStore;
0219 
0220   G4double fWeight = 1.0;
0221   G4double aScaleFactor = 1.0;
0222   G4double theLastCrossSection = 0.0;
0223   G4double mfpKinEnergy = DBL_MAX;
0224   G4int epReportLevel = 0;
0225 
0226   G4HadXSType fXSType = fHadNoIntegral;
0227 
0228 private:
0229     
0230   G4EnergyRangeManager theEnergyRangeManager;
0231   G4Nucleus targetNucleus;
0232     
0233   G4HadronicInteraction* theInteraction = nullptr;
0234   G4HadronicProcessStore* theProcessStore;
0235   const G4HadronicProcess* masterProcess = nullptr;
0236   const G4ParticleDefinition* firstParticle = nullptr;
0237   const G4ParticleDefinition* currentParticle = nullptr;
0238   const G4Material* currentMat = nullptr;
0239   const G4DynamicParticle* fDynParticle = nullptr;
0240 
0241   std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
0242   std::vector<G4TwoPeaksHadXS*>* fXSpeaks = nullptr;
0243 
0244   G4double theMFP = DBL_MAX;
0245   G4double minKinEnergy;
0246 
0247   // counters
0248   G4int nMatWarn = 0;
0249   G4int nKaonWarn = 0;
0250   G4int nICelectrons = 0;
0251   G4int matIdx = 0;
0252 
0253   // flags
0254   G4bool levelsSetByProcess = false;
0255   G4bool useIntegralXS = true;
0256   G4bool isMaster = true;
0257 
0258   G4ThreeVector unitVector;
0259 
0260   // Energy-momentum checking
0261   std::pair<G4double, G4double> epCheckLevels;
0262   std::vector<G4VLeadingParticleBiasing*> theBias;
0263 };
0264 
0265 inline G4double G4HadronicProcess::
0266 GetMicroscopicCrossSection(const G4DynamicParticle * part, 
0267                const G4Element * elm, 
0268                const G4Material* mat)
0269 { 
0270   return GetElementCrossSection(part, elm, mat);
0271 }
0272 
0273 inline G4HadronicInteraction* 
0274 G4HadronicProcess::GetHadronicInteraction() const
0275 { 
0276   return theInteraction;
0277 }
0278 
0279 inline const G4Nucleus*
0280 G4HadronicProcess::GetTargetNucleus() const
0281 {
0282   return &targetNucleus;
0283 }
0284   
0285 inline const G4Isotope* G4HadronicProcess::GetTargetIsotope()
0286 { 
0287   return targetNucleus.GetIsotope();
0288 }
0289 
0290 inline G4HadXSType 
0291 G4HadronicProcess::CrossSectionType() const
0292 { 
0293   return fXSType;
0294 }
0295 
0296 inline void 
0297 G4HadronicProcess::SetCrossSectionType(G4HadXSType val)
0298 { 
0299   fXSType = val;
0300 }
0301 
0302 inline G4double G4HadronicProcess::CrossSectionFactor() const
0303 { 
0304   return aScaleFactor;
0305 }
0306 
0307 inline void G4HadronicProcess::SetIntegral(G4bool val)
0308 { 
0309   useIntegralXS = val;
0310 }
0311 
0312 inline void G4HadronicProcess::SetEpReportLevel(G4int level)
0313 { 
0314   epReportLevel = level;
0315 }
0316 
0317 inline void 
0318 G4HadronicProcess::SetEnergyMomentumCheckLevels(G4double relativeLevel, 
0319                                                 G4double absoluteLevel)
0320 { 
0321   epCheckLevels.first = relativeLevel;
0322   epCheckLevels.second = absoluteLevel;
0323   levelsSetByProcess = true;
0324 }
0325 
0326 inline std::pair<G4double, G4double>
0327 G4HadronicProcess::GetEnergyMomentumCheckLevels() const
0328 { 
0329   return epCheckLevels;
0330 }
0331 
0332 inline G4CrossSectionDataStore*
0333 G4HadronicProcess::GetCrossSectionDataStore()
0334 {
0335   return theCrossSectionDataStore;
0336 }
0337 
0338 inline std::vector<G4TwoPeaksHadXS*>* 
0339 G4HadronicProcess::TwoPeaksXS() const
0340 { 
0341   return fXSpeaks;
0342 }
0343 
0344 inline std::vector<G4double>*
0345 G4HadronicProcess::EnergyOfCrossSectionMax() const
0346 {
0347   return theEnergyOfCrossSectionMax;
0348 }
0349 
0350 inline G4HadronicInteraction* G4HadronicProcess::
0351 ChooseHadronicInteraction(const G4HadProjectile& aHadProjectile,
0352                           G4Nucleus& aTargetNucleus,
0353                           const G4Material* aMaterial,
0354                           const G4Element* anElement)
0355 { 
0356   return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile, 
0357                                                       aTargetNucleus,
0358                                                       aMaterial,anElement);
0359 }
0360 
0361 inline G4Nucleus* G4HadronicProcess::GetTargetNucleusPointer() 
0362 { 
0363   return &targetNucleus;
0364 }
0365 
0366 inline G4double G4HadronicProcess::GetLastCrossSection() 
0367 { 
0368   return theLastCrossSection;
0369 }
0370 
0371 inline void G4HadronicProcess::DefineXSandMFP()
0372 {
0373   theLastCrossSection = aScaleFactor*
0374     theCrossSectionDataStore->GetCrossSection(fDynParticle, currentMat);
0375   theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
0376 }
0377 
0378 inline void G4HadronicProcess::ComputeXSandMFP()
0379 {
0380   theLastCrossSection = aScaleFactor*
0381     theCrossSectionDataStore->ComputeCrossSection(fDynParticle, currentMat);
0382   theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
0383 }
0384  
0385 #endif
0386