Geant4 Cross Reference |
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 // GEANT4 Class header file 29 // 30 // 31 // File name: G4VCrossSectionDataSet 32 // 33 // Author F.W. Jones, TRIUMF, 20-JAN-97 34 // 35 // Modifications: 36 // 23.01.2009 V.Ivanchenko move constructor and destructor to source 37 // 05.07.2010 V.Ivanchenko added name, min and max energy limit and 38 // corresponding access methods 39 // 12.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class 40 // 41 // 42 // Class Description 43 // This is a base class for hadronic cross section data sets. Users may 44 // derive specialized cross section classes and register them with the 45 // appropriate process, or use provided data sets. 46 // 47 // Each cross section should have unique name 48 // Minimal and maximal energy for the cross section will be used in run 49 // time before IsApplicable method is called 50 // 51 // Both the name and the energy interval will be used for documentation 52 // 53 // Class Description - End 54 55 #ifndef G4VCrossSectionDataSet_h 56 #define G4VCrossSectionDataSet_h 1 57 58 #include "globals.hh" 59 #include "G4ParticleDefinition.hh" 60 #include "G4DynamicParticle.hh" 61 #include "G4Element.hh" 62 #include <iostream> 63 64 class G4DynamicParticle; 65 class G4Isotope; 66 class G4Material; 67 class G4CrossSectionDataSetRegistry; 68 69 class G4VCrossSectionDataSet 70 { 71 public: //with description 72 73 G4VCrossSectionDataSet(const G4String& nam = ""); 74 75 virtual ~G4VCrossSectionDataSet(); 76 77 //============== Is Applicable methods =============================== 78 // The following three methods have default implementations returning 79 // "false". Derived classes should implement only needed methods. 80 81 // Element-wise cross section 82 virtual 83 G4bool IsElementApplicable(const G4DynamicParticle*, G4int Z, 84 const G4Material* mat = nullptr); 85 86 // Derived classes should implement this method if they provide isotope-wise 87 // cross sections. Default arguments G4Element and G4Material are needed to 88 // access low-energy neutron cross sections, but are not required for others. 89 virtual 90 G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A, 91 const G4Element* elm = nullptr, 92 const G4Material* mat = nullptr); 93 94 //============== GetCrossSection methods =============================== 95 96 // This is a generic method to access cross section per element 97 // This method should not be overwritten in a derived class 98 inline G4double GetCrossSection(const G4DynamicParticle*, const G4Element*, 99 const G4Material* mat = nullptr); 100 101 // This is a generic method to compute cross section per element 102 // If the DataSet is not applicable the method returns zero 103 // This method should not be overwritten in a derived class 104 G4double ComputeCrossSection(const G4DynamicParticle*, 105 const G4Element*, 106 const G4Material* mat = nullptr); 107 108 // Implement element cross section, IsApplicable does not checked. 109 // In the default implementation a sum of isotope cross sections is computed 110 virtual 111 G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, 112 const G4ParticleDefinition*, 113 const G4Element*, 114 const G4Material* mat = nullptr); 115 116 // Implement these methods for element-wise cross section 117 virtual 118 G4double GetElementCrossSection(const G4DynamicParticle*, G4int Z, 119 const G4Material* mat = nullptr); 120 121 // Derived classes should implement these methods if they provide isotope-wise 122 // cross sections. Extra arguments G4Isotope, G4Element, and G4Material are 123 // needed to access low-energy neutron cross sections, but not in other cases. 124 virtual 125 G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A, 126 const G4Isotope* iso = nullptr, 127 const G4Element* elm = nullptr, 128 const G4Material* mat = nullptr); 129 130 virtual 131 G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, 132 const G4ParticleDefinition*, G4int Z, G4int A, 133 const G4Isotope* iso = nullptr, 134 const G4Element* elm = nullptr, 135 const G4Material* mat = nullptr); 136 137 //===================================================================== 138 139 // Implement this method if needed 140 // This method is called for element-wise cross section 141 // Default implementation assumes equal cross sections for all isotopes 142 virtual const G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy, 143 G4double logE); 144 145 // Implement this method if needed 146 virtual 147 void BuildPhysicsTable(const G4ParticleDefinition&); 148 149 // Implement this method if needed 150 // Default implementation will provide a dump of the cross section 151 // in logarithmic scale in the interval of applicability 152 virtual 153 void DumpPhysicsTable(const G4ParticleDefinition&); 154 155 virtual void CrossSectionDescription(std::ostream&) const; 156 157 virtual void SetVerboseLevel(G4int value); 158 159 inline G4double GetMinKinEnergy() const; 160 161 inline void SetMinKinEnergy(G4double value); 162 163 inline G4double GetMaxKinEnergy() const; 164 165 inline void SetMaxKinEnergy(G4double value); 166 167 inline bool ForAllAtomsAndEnergies() const; 168 169 inline void SetForAllAtomsAndEnergies(G4bool val); 170 171 inline const G4String& GetName() const; 172 173 inline void SetName(const G4String& nam); 174 175 G4VCrossSectionDataSet & operator= 176 (const G4VCrossSectionDataSet &right) = delete; 177 G4VCrossSectionDataSet(const G4VCrossSectionDataSet&) = delete; 178 179 protected: 180 181 G4int verboseLevel{0}; 182 183 G4String name; 184 185 private: 186 187 G4CrossSectionDataSetRegistry* registry; 188 189 G4double minKinEnergy{0.0}; 190 G4double maxKinEnergy; 191 192 G4bool isForAllAtomsAndEnergies{false}; 193 194 }; 195 196 inline G4double 197 G4VCrossSectionDataSet::GetCrossSection(const G4DynamicParticle* dp, 198 const G4Element* elm, 199 const G4Material* mat) 200 { 201 return ComputeCrossSection(dp, elm, mat); 202 } 203 204 inline void G4VCrossSectionDataSet::SetVerboseLevel(G4int value) 205 { 206 verboseLevel = value; 207 } 208 209 inline void G4VCrossSectionDataSet::SetMinKinEnergy(G4double value) 210 { 211 minKinEnergy = value; 212 } 213 214 inline G4double G4VCrossSectionDataSet::GetMinKinEnergy() const 215 { 216 return minKinEnergy; 217 } 218 219 inline void G4VCrossSectionDataSet::SetMaxKinEnergy(G4double value) 220 { 221 maxKinEnergy = value; 222 } 223 224 inline G4double G4VCrossSectionDataSet::GetMaxKinEnergy() const 225 { 226 return maxKinEnergy; 227 } 228 229 inline const G4String& G4VCrossSectionDataSet::GetName() const 230 { 231 return name; 232 } 233 234 inline bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies() const 235 { 236 return isForAllAtomsAndEnergies; 237 } 238 239 inline void G4VCrossSectionDataSet::SetForAllAtomsAndEnergies(G4bool val) 240 { 241 isForAllAtomsAndEnergies = val; 242 } 243 244 inline void G4VCrossSectionDataSet::SetName(const G4String& nam) 245 { 246 name = nam; 247 } 248 249 #endif 250