Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VCrossSectionDataSet.hh
이 파일의 문서화 페이지로 가기
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 // $Id: $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4VCrossSectionDataSet
34 //
35 // Author F.W. Jones, TRIUMF, 20-JAN-97
36 //
37 // Modifications:
38 // 23.01.2009 V.Ivanchenko move constructor and destructor to source
39 // 05.07.2010 V.Ivanchenko added name, min and max energy limit and
40 // corresponding access methods
41 // 12.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class
42 //
43 //
44 // Class Description
45 // This is a base class for hadronic cross section data sets. Users may
46 // derive specialized cross section classes and register them with the
47 // appropriate process, or use provided data sets.
48 //
49 // Each cross section should have unique name
50 // Minimal and maximal energy for the cross section will be used in run
51 // time before IsApplicable method is called
52 //
53 // Both the name and the energy interval will be used for documentation
54 //
55 // Class Description - End
56 
57 #ifndef G4VCrossSectionDataSet_h
58 #define G4VCrossSectionDataSet_h 1
59 
60 #include "globals.hh"
61 #include "G4ParticleDefinition.hh"
62 #include "G4Element.hh"
63 #include "G4HadTmpUtil.hh"
64 #include <iostream>
65 
66 class G4DynamicParticle;
67 class G4Isotope;
68 class G4Material;
70 
72 {
73 public: //with description
74 
75  G4VCrossSectionDataSet(const G4String& nam = "");
76 
77  virtual ~G4VCrossSectionDataSet();
78 
79  //============== Is Applicable methods ===============================
80  // The following three methods have default implementations returning
81  // "false". Derived classes should implement only needed methods.
82 
83  // Element-wise cross section
84  virtual
86  const G4Material* mat = 0);
87 
88  // Derived classes should implement this method if they provide isotope-wise
89  // cross sections. Default arguments G4Element and G4Material are needed to
90  // access low-energy neutron cross sections, but are not required for others.
91  virtual
93  const G4Element* elm = 0,
94  const G4Material* mat = 0);
95 
96  //============== GetCrossSection methods ===============================
97 
98  // This is a generic method to access cross section per element
99  // This method should not be overwritten in a derived class
100  inline G4double GetCrossSection(const G4DynamicParticle*, const G4Element*,
101  const G4Material* mat = nullptr);
102 
103  // This is a generic method to compute cross section per element
104  // If the DataSet is not applicable the method returns zero
105  // This method should not be overwritten in a derived class
107  const G4Element*,
108  const G4Material* mat = nullptr);
109 
110  // The following two methods have default implementations which throw
111  // G4HadronicException. Derived classes should implement only needed
112  // methods, which are assumed to be called at run time.
113 
114  // Implement this method for element-wise cross section
115  virtual
117  const G4Material* mat = nullptr);
118 
119  // Derived classes should implement this method if they provide isotope-wise
120  // cross sections. Default arguments G4Element and G4Material are needed to
121  // access low-energy neutron cross sections, but are not required for others.
122  virtual
124  const G4Isotope* iso = nullptr,
125  const G4Element* elm = nullptr,
126  const G4Material* mat = nullptr);
127 
128  //=====================================================================
129 
130  // Implement this method if needed
131  // This method is called for element-wise cross section
132  // Default implementation assumes equal cross sections for all isotopes
133  virtual const G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy);
134 
135  // Implement this method if needed
136  virtual
138 
139  // Implement this method if needed
140  // Default implementation will provide a dump of the cross section
141  // in logarithmic scale in the interval of applicability
142  virtual
144 
145  virtual void CrossSectionDescription(std::ostream&) const;
146 
147 public: // Without Description
148 
149  virtual G4int GetVerboseLevel() const;
150 
151  virtual void SetVerboseLevel(G4int value);
152 
153  inline G4double GetMinKinEnergy() const;
154 
155  inline void SetMinKinEnergy(G4double value);
156 
157  inline G4double GetMaxKinEnergy() const;
158 
159  inline void SetMaxKinEnergy(G4double value);
160 
161  inline bool ForAllAtomsAndEnergies() const;
162 
163  inline void SetForAllAtomsAndEnergies(G4bool val);
164 
165  inline const G4String& GetName() const;
166 
167 protected:
168 
169  inline void SetName(const G4String&);
170 
172 
173 private:
174 
177 
179 
182 
184 
186 };
187 
188 inline G4double
190  const G4Element* elm,
191  const G4Material* mat)
192 {
193  return ComputeCrossSection(dp, elm, mat);
194 }
195 
196 
198 {
199  return verboseLevel;
200 }
201 
203 {
205 }
206 
208 {
210 }
211 
213 {
214  return minKinEnergy;
215 }
216 
218 {
220 }
221 
223 {
224  return maxKinEnergy;
225 }
226 
228 {
229  return name;
230 }
231 
233 {
235 }
236 
238 {
240 }
241 
243 {
244  name = nam;
245 }
246 
247 #endif
const XML_Char * name
Definition: expat.h:151
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
void SetMinKinEnergy(G4double value)
G4CrossSectionDataSetRegistry * registry
void SetMaxKinEnergy(G4double value)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
virtual void SetVerboseLevel(G4int value)
const G4String & GetName() const
Float_t Z
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
virtual const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
virtual G4int GetVerboseLevel() const
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4VCrossSectionDataSet(const G4String &nam="")
double A(double temperature)
void SetForAllAtomsAndEnergies(G4bool val)
Float_t mat
int G4int
Definition: G4Types.hh:78
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
void SetName(const G4String &)
virtual void CrossSectionDescription(std::ostream &) const
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
G4VCrossSectionDataSet & operator=(const G4VCrossSectionDataSet &right)