Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/include/G4CrossSectionDataStore.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 //
 27 // -------------------------------------------------------------------
 28 // File name:     G4CrossSectionDataStore
 29 //
 30 // Modifications:
 31 // 23.01.2009 V.Ivanchenko move constructor and destructor to source,
 32 //                         use STL vector instead of C-array
 33 //
 34 // August 2011  Re-designed
 35 //              by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
 36 
 37 // Class Description
 38 // This is the class to which cross section data sets may be registered. 
 39 // An instance of it is contained in each hadronic process, allowing
 40 // the use of the AddDataSet() method to tailor the cross sections to
 41 // your application.
 42 // Class Description - End
 43 
 44 #ifndef G4CrossSectionDataStore_h
 45 #define G4CrossSectionDataStore_h 1
 46 
 47 #include "globals.hh"
 48 #include "G4VCrossSectionDataSet.hh"
 49 #include "G4DynamicParticle.hh"
 50 #include "G4PhysicsVector.hh"
 51 #include <vector>
 52 #include <iostream>
 53 
 54 class G4Nucleus;
 55 class G4ParticleDefinition;
 56 class G4Isotope;
 57 class G4Element;
 58 class G4Material;
 59 class G4NistManager;
 60 
 61 class G4CrossSectionDataStore
 62 {
 63 public:
 64 
 65   G4CrossSectionDataStore();
 66 
 67   ~G4CrossSectionDataStore() = default;
 68 
 69   // Run time cross section per unit volume
 70   inline G4double GetCrossSection(const G4DynamicParticle*, const G4Material*);
 71   G4double ComputeCrossSection(const G4DynamicParticle*, const G4Material*);
 72 
 73   // Cross section per element
 74   G4double GetCrossSection(const G4DynamicParticle*, 
 75          const G4Element*, const G4Material*);
 76 
 77   // Cross section per isotope 
 78   G4double GetCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
 79                            const G4Isotope*,
 80          const G4Element*, const G4Material*);
 81 
 82   // Sample Z and A of a target nucleus and upload into G4Nucleus
 83   const G4Element* SampleZandA(const G4DynamicParticle*, const G4Material*,
 84              G4Nucleus& target);
 85 
 86   // Initialisation before run
 87   void BuildPhysicsTable(const G4ParticleDefinition&);
 88 
 89   // Dump store to G4cout
 90   void DumpPhysicsTable(const G4ParticleDefinition&);
 91 
 92   // Dump store as html
 93   void DumpHtml(const G4ParticleDefinition&, std::ofstream&) const;
 94   void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs,
 95                              const G4String&, const G4String&) const;
 96   
 97   void AddDataSet(G4VCrossSectionDataSet*);
 98   void AddDataSet(G4VCrossSectionDataSet*, std::size_t);
 99   inline const std::vector<G4VCrossSectionDataSet*>& GetDataSetList() const;
100 
101   inline void SetVerboseLevel(G4int value);
102 
103   // may be used by special processes
104   inline void SetForcedElement(const G4Element*);
105 
106   G4CrossSectionDataStore & operator=
107   (const G4CrossSectionDataStore &right) = delete;
108   G4CrossSectionDataStore(const G4CrossSectionDataStore&) = delete;
109 
110 private:
111 
112   G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
113             const G4Isotope*, const G4Element*,
114                               const G4Material*, const G4int index);
115 
116   G4String HtmlFileName(const G4String & in) const;
117 
118   G4NistManager* nist;
119   const G4Material* currentMaterial = nullptr;
120   const G4ParticleDefinition* matParticle = nullptr;
121   const G4Element* forcedElement = nullptr;
122   G4double matKinEnergy = 0.0;
123   G4double matCrossSection = 0.0;
124 
125   G4int nDataSetList = 0;
126   G4int verboseLevel = 1;
127 
128   std::vector<G4VCrossSectionDataSet*> dataSetList;
129   std::vector<G4double> xsecelm;
130   std::vector<G4double> xseciso;
131 };
132 
133 inline void G4CrossSectionDataStore::SetVerboseLevel(G4int value)
134 {
135   verboseLevel = value;
136 }
137 
138 inline void G4CrossSectionDataStore::SetForcedElement(const G4Element* ptr)
139 {
140   forcedElement = ptr;
141 }
142 
143 inline const std::vector<G4VCrossSectionDataSet*>&
144 G4CrossSectionDataStore::GetDataSetList() const
145 {
146   return dataSetList;
147 }
148 
149 inline G4double 
150 G4CrossSectionDataStore::GetCrossSection(const G4DynamicParticle* dp,
151                                          const G4Material* mat)
152 {
153   if(dp->GetKineticEnergy() != matKinEnergy || mat != currentMaterial ||
154      dp->GetDefinition() != matParticle) {
155     ComputeCrossSection(dp, mat);
156   }
157   return matCrossSection;
158 }
159 
160 #endif
161