Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PhysicsVector.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// G4PhysicsVector
27//
28// Class description:
29//
30// A physics vector which has values of energy-loss, cross-section,
31// and other physics values of a particle in matter in a given
32// range of energy, momentum, etc.
33// This class serves as the base class for a vector having various
34// energy scale, for example like 'log', 'linear', 'free', etc.
35
36// Authors:
37// - 02 Dec. 1995, G.Cosmo: Structure created based on object model
38// - 03 Mar. 1996, K.Amako: Implemented the 1st version
39// Revisions:
40// - 11 Nov. 2000, H.Kurashige: Use STL vector for dataVector and binVector
41// --------------------------------------------------------------------
42#ifndef G4PhysicsVector_hh
43#define G4PhysicsVector_hh 1
44
45#include <fstream>
46#include <iostream>
47#include <vector>
48
49#include "G4Log.hh"
51#include "G4ios.hh"
52#include "globals.hh"
53
55{
56public:
57 // Default constructor - vector will be filled via Retrieve() method
58 // Free vector may be filled via InsertValue(..) method
59 explicit G4PhysicsVector(G4bool spline = false);
60
61 // Copy constructor and assignment operator
64
65 // not used operators
68 G4bool operator==(const G4PhysicsVector& right) const = delete;
69 G4bool operator!=(const G4PhysicsVector& right) const = delete;
70
71 virtual ~G4PhysicsVector() = default;
72
73 // Get the cross-section/energy-loss value corresponding to the
74 // given energy. An appropriate interpolation is used to calculate
75 // the value. Consumer code gets changed index and may reuse it
76 // for the next call to save CPU for bin location.
77 inline G4double Value(const G4double energy, std::size_t& lastidx) const;
78
79 // Get the cross-section/energy-loss value corresponding to the
80 // given energy. An appropriate interpolation is used to calculate
81 // the value. This method should be used if bin location cannot be
82 // kept in the user code.
83 inline G4double Value(const G4double energy) const;
84
85 // Obsolete method to get value, 'isOutRange' is not used anymore.
86 // This method is kept for the compatibility reason
87 inline G4double GetValue(const G4double energy, G4bool& isOutRange) const;
88
89 // Same as the Value() method above but specialised for log-vector type.
90 // Note, unlike the general Value() method above, this method will work
91 // properly only for G4PhysicsLogVector.
92 inline G4double LogVectorValue(const G4double energy,
93 const G4double theLogEnergy) const;
94
95 // Returns the value for the specified index of the dataVector
96 // The boundary check will not be done
97 inline G4double operator[](const std::size_t index) const;
98 inline G4double operator()(const std::size_t index) const;
99
100 // Put data into the vector at 'index' position.
101 // Take note that the 'index' starts from '0'.
102 // It is assumed that energies are already filled.
103 inline void PutValue(const std::size_t index, const G4double value);
104
105 // Returns the value in the energy specified by 'index'
106 // of the energy vector. The boundary check will not be done.
107 // Use this when compute cross-section, dEdx, or other value
108 // before filling the vector by PutValue().
109 inline G4double Energy(const std::size_t index) const;
110 inline G4double GetLowEdgeEnergy(const std::size_t index) const;
111
112 // Returns the energy of the first and the last point of the vector.
113 inline G4double GetMinEnergy() const;
114 inline G4double GetMaxEnergy() const;
115
116 // Returns the data of the first and the last point of the vector.
117 // If the vector is empty returns zeros.
118 inline G4double GetMinValue() const;
119 inline G4double GetMaxValue() const;
120
121 // Get the total length of the vector
122 inline std::size_t GetVectorLength() const;
123
124 // Computes the lower index the energy bin in case of log-vector i.e.
125 // in case of vectors with equal bin widths on log-scale
126 // Note, that no check on the boundary is performed
127 inline std::size_t ComputeLogVectorBin(const G4double logenergy) const;
128
129 // Get physics vector type.
131
132 // True if using spline interpolation.
133 inline G4bool GetSpline() const;
134
135 // Define verbosity level.
136 inline void SetVerboseLevel(G4int value);
137
138 // Find energy using linear interpolation for vector
139 // filled by cumulative probability function.
140 // Assuming that vector is already filled.
141 inline G4double FindLinearEnergy(const G4double rand) const;
142
143 // Find low edge index of a bin for given energy.
144 // Min value 0, max value idxmax.
145 std::size_t FindBin(const G4double energy, std::size_t idx) const;
146
147 // Scale all values of the vector by factorV, energies by vectorE.
148 // AFter this method FillSecondDerivatives(...) should be called.
149 // This method may be applied for example after retrieving a vector
150 // from an external file to convert values into Geant4 units.
151 void ScaleVector(const G4double factorE, const G4double factorV);
152
153 // This method should be called when the vector is fully filled
154 // There are 3 types of second derivative computations:
155 // fSplineSimple - 2d derivative continues
156 // fSplineBase - 3d derivative continues (the default)
157 // fSplineFixedEdges - 3d derivatives continues, 1st and last
158 // derivatives are fixed
159 void FillSecondDerivatives(const G4SplineType = G4SplineType::Base,
160 const G4double dir1 = 0.0,
161 const G4double dir2 = 0.0);
162
163 // This method can be applied if both energy and data values
164 // grow monotonically, for example, if in this vector a
165 // cumulative probability density function is stored.
166 G4double GetEnergy(const G4double value) const;
167
168 // To store/retrieve persistent data to/from file streams.
169 G4bool Store(std::ofstream& fOut, G4bool ascii = false) const;
170 G4bool Retrieve(std::ifstream& fIn, G4bool ascii = false);
171
172 // Print vector
173 friend std::ostream& operator<<(std::ostream&, const G4PhysicsVector&);
174 void DumpValues(G4double unitE = 1.0, G4double unitV = 1.0) const;
175
176protected:
177
178 // The default implements a free vector initialisation.
179 virtual void Initialise();
180
181 void PrintPutValueError(std::size_t index, G4double value,
182 const G4String& text);
183
184private:
185
186 void ComputeSecDerivative0();
187 void ComputeSecDerivative1();
188 void ComputeSecDerivative2(const G4double firstPointDerivative,
189 const G4double endPointDerivative);
190 // Internal methods for computing of spline coeffitients
191
192 // Linear or spline interpolation.
193 inline G4double Interpolation(const std::size_t idx,
194 const G4double energy) const;
195
196 // Assuming (edgeMin <= energy <= edgeMax).
197 inline std::size_t GetBin(const G4double energy) const;
198
199protected:
200
201 G4double edgeMin = 0.0; // Energy of first point
202 G4double edgeMax = 0.0; // Energy of the last point
203
204 G4double invdBin = 0.0; // 1/Bin width for linear and log vectors
205 G4double logemin = 0.0; // used only for log vector
206
208 std::size_t idxmax = 0;
209 std::size_t numberOfNodes = 0;
210
212 // The type of PhysicsVector (enumerator)
213
214 std::vector<G4double> binVector; // energy
215 std::vector<G4double> dataVector; // crossection/energyloss
216 std::vector<G4double> secDerivative; // second derivatives
217
218private:
219
220 G4bool useSpline = false;
221};
222
223#include "G4PhysicsVector.icc"
224
225#endif
G4PhysicsVectorType
@ T_G4PhysicsFreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetVerboseLevel(G4int value)
G4bool operator==(const G4PhysicsVector &right) const =delete
G4double Value(const G4double energy) const
std::size_t idxmax
virtual ~G4PhysicsVector()=default
G4PhysicsVectorType type
G4double GetEnergy(const G4double value) const
G4PhysicsVector & operator=(const G4PhysicsVector &&)=delete
G4double GetMinValue() const
G4PhysicsVector(const G4PhysicsVector &&)=delete
G4double GetLowEdgeEnergy(const std::size_t index) const
void PutValue(const std::size_t index, const G4double value)
void PrintPutValueError(std::size_t index, G4double value, const G4String &text)
G4double GetMaxEnergy() const
G4double operator[](const std::size_t index) const
std::size_t ComputeLogVectorBin(const G4double logenergy) const
void ScaleVector(const G4double factorE, const G4double factorV)
G4double GetValue(const G4double energy, G4bool &isOutRange) const
G4bool Store(std::ofstream &fOut, G4bool ascii=false) const
G4PhysicsVector(const G4PhysicsVector &)=default
G4double GetMaxValue() const
G4double Energy(const std::size_t index) const
std::size_t numberOfNodes
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const
G4PhysicsVectorType GetType() const
std::vector< G4double > secDerivative
G4double operator()(const std::size_t index) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
std::vector< G4double > dataVector
std::vector< G4double > binVector
G4double GetMinEnergy() const
virtual void Initialise()
std::size_t GetVectorLength() const
G4bool operator!=(const G4PhysicsVector &right) const =delete
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
friend std::ostream & operator<<(std::ostream &, const G4PhysicsVector &)
G4double FindLinearEnergy(const G4double rand) const
G4PhysicsVector & operator=(const G4PhysicsVector &)=default
G4bool GetSpline() const
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
std::size_t FindBin(const G4double energy, std::size_t idx) const