00001 // 00002 // ******************************************************************** 00003 // * License and Disclaimer * 00004 // * * 00005 // * The Geant4 software is copyright of the Copyright Holders of * 00006 // * the Geant4 Collaboration. It is provided under the terms and * 00007 // * conditions of the Geant4 Software License, included in the file * 00008 // * LICENSE and available at http://cern.ch/geant4/license . These * 00009 // * include a list of copyright holders. * 00010 // * * 00011 // * Neither the authors of this software system, nor their employing * 00012 // * institutes,nor the agencies providing financial support for this * 00013 // * work make any representation or warranty, express or implied, * 00014 // * regarding this software system or assume any liability for its * 00015 // * use. Please see the license in the file LICENSE and URL above * 00016 // * for the full disclaimer and the limitation of liability. * 00017 // * * 00018 // * This code implementation is the result of the scientific and * 00019 // * technical work of the GEANT4 collaboration. * 00020 // * By using, copying, modifying or distributing the software (or * 00021 // * any work based on the software) you agree to acknowledge its * 00022 // * use in resulting scientific publications, and indicate your * 00023 // * acceptance of all terms of the Geant4 Software license. * 00024 // ******************************************************************** 00025 // 00026 // $Id: $ 00027 // GEANT4 tag $Name: not supported by cvs2svn $ 00028 // 00029 // ------------------------------------------------------------------- 00030 // 00031 // GEANT4 Class header file 00032 // 00033 // 00034 // File name: G4VCrossSectionDataSet 00035 // 00036 // Author F.W. Jones, TRIUMF, 20-JAN-97 00037 // 00038 // Modifications: 00039 // 23.01.2009 V.Ivanchenko move constructor and destructor to source 00040 // 05.07.2010 V.Ivanchenko added name, min and max energy limit and 00041 // corresponding access methods 00042 // 12.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class 00043 // 00044 // 00045 // Class Description 00046 // This is a base class for hadronic cross section data sets. Users may 00047 // derive specialized cross section classes and register them with the 00048 // appropriate process, or use provided data sets. 00049 // 00050 // Each cross section should have unique name 00051 // Minimal and maximal energy for the cross section will be used in run 00052 // time before IsApplicable method is called 00053 // 00054 // Both the name and the energy interval will be used for documentation 00055 // 00056 // Class Description - End 00057 00058 #ifndef G4VCrossSectionDataSet_h 00059 #define G4VCrossSectionDataSet_h 1 00060 00061 #include "globals.hh" 00062 #include "G4ParticleDefinition.hh" 00063 #include "G4Element.hh" 00064 #include "G4HadTmpUtil.hh" 00065 #include <iostream> 00066 00067 class G4DynamicParticle; 00068 class G4Isotope; 00069 class G4Material; 00070 00071 class G4VCrossSectionDataSet 00072 { 00073 public: //with description 00074 00075 G4VCrossSectionDataSet(const G4String& nam = ""); 00076 00077 virtual ~G4VCrossSectionDataSet(); 00078 00079 //============== Is Applicable methods =============================== 00080 // The following three methods have default implementations returning 00081 // "false". Derived classes should implement only needed methods. 00082 00083 // Element-wise cross section 00084 virtual 00085 G4bool IsElementApplicable(const G4DynamicParticle*, G4int Z, 00086 const G4Material* mat = 0); 00087 00088 // Derived classes should implement this method if they provide isotope-wise 00089 // cross sections. Default arguments G4Element and G4Material are needed to 00090 // access low-energy neutron cross sections, but are not required for others. 00091 virtual 00092 G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A, 00093 const G4Element* elm = 0, 00094 const G4Material* mat = 0); 00095 00096 //============== GetCrossSection methods =============================== 00097 00098 // This is a generic method to access cross section per element 00099 // This method should not be overwritten in a derived class 00100 inline G4double GetCrossSection(const G4DynamicParticle*, const G4Element*, 00101 const G4Material* mat = 0); 00102 00103 // This is a generic method to compute cross section per element 00104 // If the DataSet is not applicable the method returns zero 00105 // This method should not be overwritten in a derived class 00106 G4double ComputeCrossSection(const G4DynamicParticle*, 00107 const G4Element*, 00108 const G4Material* mat = 0); 00109 00110 // The following two methods have default implementations which throw 00111 // G4HadronicException. Derived classes should implement only needed 00112 // methods, which are assumed to be called at run time. 00113 00114 // Implement this method for element-wise cross section 00115 virtual 00116 G4double GetElementCrossSection(const G4DynamicParticle*, G4int Z, 00117 const G4Material* mat = 0); 00118 00119 // Derived classes should implement this method if they provide isotope-wise 00120 // cross sections. Default arguments G4Element and G4Material are needed to 00121 // access low-energy neutron cross sections, but are not required for others. 00122 virtual 00123 G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A, 00124 const G4Isotope* iso = 0, 00125 const G4Element* elm = 0, 00126 const G4Material* mat = 0); 00127 00128 //===================================================================== 00129 00130 // Implement this method if needed 00131 // This method is called for element-wise cross section 00132 // Default implementation assumes equal cross sections for all isotopes 00133 virtual G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy); 00134 00135 // Implement this method if needed 00136 virtual 00137 void BuildPhysicsTable(const G4ParticleDefinition&); 00138 00139 // Implement this method if needed 00140 // Default implementation will provide a dump of the cross section 00141 // in logarithmic scale in the interval of applicability 00142 virtual 00143 void DumpPhysicsTable(const G4ParticleDefinition&); 00144 00145 virtual void CrossSectionDescription(std::ostream&) const; 00146 00147 public: // Without Description 00148 00149 inline void SetVerboseLevel(G4int value); 00150 00151 inline G4double GetMinKinEnergy() const; 00152 00153 inline void SetMinKinEnergy(G4double value); 00154 00155 inline G4double GetMaxKinEnergy() const; 00156 00157 inline void SetMaxKinEnergy(G4double value); 00158 00159 inline const G4String& GetName() const; 00160 00161 protected: 00162 00163 inline void SetName(const G4String&); 00164 00165 G4int verboseLevel; 00166 00167 private: 00168 00169 G4VCrossSectionDataSet & operator=(const G4VCrossSectionDataSet &right); 00170 G4VCrossSectionDataSet(const G4VCrossSectionDataSet&); 00171 00172 G4double minKinEnergy; 00173 G4double maxKinEnergy; 00174 00175 G4String name; 00176 }; 00177 00178 inline G4double 00179 G4VCrossSectionDataSet::GetCrossSection(const G4DynamicParticle* dp, 00180 const G4Element* elm, 00181 const G4Material* mat) 00182 { 00183 return ComputeCrossSection(dp, elm, mat); 00184 } 00185 00186 inline void G4VCrossSectionDataSet::SetVerboseLevel(G4int value) 00187 { 00188 verboseLevel = value; 00189 } 00190 00191 inline void G4VCrossSectionDataSet::SetMinKinEnergy(G4double value) 00192 { 00193 minKinEnergy = value; 00194 } 00195 00196 inline G4double G4VCrossSectionDataSet::GetMinKinEnergy() const 00197 { 00198 return minKinEnergy; 00199 } 00200 00201 inline void G4VCrossSectionDataSet::SetMaxKinEnergy(G4double value) 00202 { 00203 maxKinEnergy = value; 00204 } 00205 00206 inline G4double G4VCrossSectionDataSet::GetMaxKinEnergy() const 00207 { 00208 return maxKinEnergy; 00209 } 00210 00211 inline const G4String& G4VCrossSectionDataSet::GetName() const 00212 { 00213 return name; 00214 } 00215 00216 inline void G4VCrossSectionDataSet::SetName(const G4String& nam) 00217 { 00218 name = nam; 00219 } 00220 00221 #endif