Geant4_10
G4VCrossSectionDataSet.hh
Go to the documentation of this file.
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;
69 
71 {
72 public: //with description
73 
74  G4VCrossSectionDataSet(const G4String& nam = "");
75 
76  virtual ~G4VCrossSectionDataSet();
77 
78  //============== Is Applicable methods ===============================
79  // The following three methods have default implementations returning
80  // "false". Derived classes should implement only needed methods.
81 
82  // Element-wise cross section
83  virtual
85  const G4Material* mat = 0);
86 
87  // Derived classes should implement this method if they provide isotope-wise
88  // cross sections. Default arguments G4Element and G4Material are needed to
89  // access low-energy neutron cross sections, but are not required for others.
90  virtual
92  const G4Element* elm = 0,
93  const G4Material* mat = 0);
94 
95  //============== GetCrossSection methods ===============================
96 
97  // This is a generic method to access cross section per element
98  // This method should not be overwritten in a derived class
100  const G4Material* mat = 0);
101 
102  // This is a generic method to compute cross section per element
103  // If the DataSet is not applicable the method returns zero
104  // This method should not be overwritten in a derived class
106  const G4Element*,
107  const G4Material* mat = 0);
108 
109  // The following two methods have default implementations which throw
110  // G4HadronicException. Derived classes should implement only needed
111  // methods, which are assumed to be called at run time.
112 
113  // Implement this method for element-wise cross section
114  virtual
116  const G4Material* mat = 0);
117 
118  // Derived classes should implement this method if they provide isotope-wise
119  // cross sections. Default arguments G4Element and G4Material are needed to
120  // access low-energy neutron cross sections, but are not required for others.
121  virtual
123  const G4Isotope* iso = 0,
124  const G4Element* elm = 0,
125  const G4Material* mat = 0);
126 
127  //=====================================================================
128 
129  // Implement this method if needed
130  // This method is called for element-wise cross section
131  // Default implementation assumes equal cross sections for all isotopes
132  virtual G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy);
133 
134  // Implement this method if needed
135  virtual
137 
138  // Implement this method if needed
139  // Default implementation will provide a dump of the cross section
140  // in logarithmic scale in the interval of applicability
141  virtual
143 
144  virtual void CrossSectionDescription(std::ostream&) const;
145 
146 public: // Without Description
147 
148  virtual G4int GetVerboseLevel() const;
149 
150  virtual void SetVerboseLevel(G4int value);
151 
152  inline G4double GetMinKinEnergy() const;
153 
154  inline void SetMinKinEnergy(G4double value);
155 
156  inline G4double GetMaxKinEnergy() const;
157 
158  inline void SetMaxKinEnergy(G4double value);
159 
160  inline const G4String& GetName() const;
161 
162 protected:
163 
164  inline void SetName(const G4String&);
165 
167 
168 private:
169 
172 
173  G4double minKinEnergy;
174  G4double maxKinEnergy;
175 
176  G4String name;
177 };
178 
179 inline G4double
181  const G4Element* elm,
182  const G4Material* mat)
183 {
184  return ComputeCrossSection(dp, elm, mat);
185 }
186 
187 
189 {
190  return verboseLevel;
191 }
192 
194 {
196 }
197 
199 {
200  minKinEnergy = value;
201 }
202 
204 {
205  return minKinEnergy;
206 }
207 
209 {
210  maxKinEnergy = value;
211 }
212 
214 {
215  return maxKinEnergy;
216 }
217 
219 {
220  return name;
221 }
222 
224 {
225  name = nam;
226 }
227 
228 #endif
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
virtual void CrossSectionDescription(std::ostream &) const
G4VCrossSectionDataSet(const G4String &nam="")
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const XML_Char * name
Definition: expat.h:151
const G4String & GetName() const
int G4int
Definition: G4Types.hh:78
Float_t mat
Definition: plot.C:40
void SetName(const G4String &)
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
Float_t Z
Definition: plot.C:39
void SetMinKinEnergy(G4double value)
bool G4bool
Definition: G4Types.hh:79
virtual G4int GetVerboseLevel() const
void SetMaxKinEnergy(G4double value)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
virtual void SetVerboseLevel(G4int value)
const XML_Char int const XML_Char * value
Definition: expat.h:331
double G4double
Definition: G4Types.hh:76
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
virtual G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy)