Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/particle_hp/include/G4CrossSectionHP.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 // V. Ivanchenko September 2023 
 27 //               
 28 // G4CrossSectionHP is a generic class implementing 
 29 // cross sections for neutrons, protons and light ions
 30 // It is an alternative to code developed by J.P. Wellisch & T.Koi
 31 //
 32 
 33 #ifndef G4CrossSectionHP_h
 34 #define G4CrossSectionHP_h 1
 35 
 36 #include "G4VCrossSectionDataSet.hh"
 37 #include "globals.hh"
 38 #include "G4ElementData.hh"
 39 #include "G4PhysicsVector.hh"
 40 #include "G4LorentzVector.hh"
 41 #include "G4ThreeVector.hh"
 42 #include <vector>
 43 
 44 class G4DynamicParticle;
 45 class G4ParticleDefinition;
 46 class G4ParticleHPManager;
 47 class G4Element;
 48 class G4Material;
 49 
 50 class G4CrossSectionHP : public G4VCrossSectionDataSet
 51 {
 52   public:
 53     explicit G4CrossSectionHP(const G4ParticleDefinition*,
 54                               const G4String& nameData,
 55                               const G4String& nameDir, G4double emaxHP,
 56                               G4int zmin, G4int zmax);
 57 
 58     ~G4CrossSectionHP() override = default;
 59 
 60     G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A,
 61                            const G4Element*, const G4Material*) override;
 62 
 63     G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
 64                                 const G4Isotope*, const G4Element*,
 65                                 const G4Material*) override;
 66 
 67     G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge,
 68                                     const G4ParticleDefinition*,
 69                                     G4int Z, G4int A,
 70                                     const G4Isotope* iso,
 71                                     const G4Element* elm,
 72                                     const G4Material* mat) override;
 73 
 74     const G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy,
 75                                    G4double logE) override;
 76 
 77     void BuildPhysicsTable(const G4ParticleDefinition&) override;
 78 
 79     void DumpPhysicsTable(const G4ParticleDefinition&) override;
 80 
 81     G4CrossSectionHP & operator=(const G4CrossSectionHP &right) = delete;
 82     G4CrossSectionHP(const G4CrossSectionHP&) = delete;
 83 
 84   protected:
 85 
 86     inline G4double GetMaxHPEnergy() const;
 87 
 88     inline void SetBinSearch(G4int n);
 89 
 90   private:
 91 
 92     void Initialise(const G4int Z);
 93 
 94     void PrepareCache(const G4Material*);
 95 
 96     G4double IsoCrossSection(const G4double kinE, const G4double loge,
 97            const G4int Z, const G4int A,
 98                              const G4double temperature);
 99 
100     inline G4double GetCrossSection(const G4int Z, const G4int A);
101 
102     inline G4bool CheckCache(const G4int Z);
103 
104     const G4ParticleDefinition* fParticle;
105     const G4ParticleDefinition* fNeutron;
106     G4ParticleHPManager* fManagerHP;
107 
108     const G4double emax;
109     const G4double emaxT;
110     const G4double elimit;
111     const G4double logElimit;
112     G4LorentzVector fLV;
113     G4ThreeVector fBoost;
114 
115     const G4int minZ;
116     const G4int maxZ;
117 
118     G4int binSearch{2};
119     std::size_t index{0};
120     G4bool fPrinted{false};
121 
122     // cache
123     const G4Material* fCurrentMat{nullptr};
124     std::vector<std::pair<G4int, G4int> > fZA;
125     std::vector<G4double> fIsoXS;
126     std::vector<G4double> fTemp;
127 
128     const G4String fDataName;
129     const G4String fDataDirectory;
130     G4ElementData* fData{nullptr};
131 };
132 
133 inline G4double
134 G4CrossSectionHP::GetCrossSection(const G4int Z, const G4int A)
135 {
136   G4double res = 0.0;
137   for (std::size_t i=0; i<fZA.size(); ++i) {
138     if (Z == fZA[i].first && A == fZA[i].second) {
139       res = fIsoXS[i];
140       break;
141     }
142   }
143   return res;
144 }
145 
146 inline G4bool G4CrossSectionHP::CheckCache(const G4int Z)
147 {
148   G4bool res = false;
149   for (auto const & p : fZA) {
150     if (Z == p.first) {
151       res = true;
152       break;
153     }
154   }
155   return res;
156 }
157 
158 inline G4double G4CrossSectionHP::GetMaxHPEnergy() const
159 { 
160   return emax;
161 }
162 
163 inline void G4CrossSectionHP::SetBinSearch(G4int n)
164 {
165   if (n > 0) { binSearch = n; }
166 }
167 
168 #endif
169