Geant4-11
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
G4Physics2DVector Class Reference

#include <G4Physics2DVector.hh>

Public Member Functions

std::size_t FindBinLocationX (const G4double x, const std::size_t lastidx) const
 
std::size_t FindBinLocationY (const G4double y, const std::size_t lastidy) const
 
G4double FindLinearX (G4double rand, G4double y) const
 
G4double FindLinearX (G4double rand, G4double y, std::size_t &lastidy) const
 
 G4Physics2DVector ()
 
 G4Physics2DVector (const G4Physics2DVector &)
 
 G4Physics2DVector (std::size_t nx, std::size_t ny)
 
std::size_t GetLengthX () const
 
std::size_t GetLengthY () const
 
G4PhysicsVectorType GetType () const
 
G4double GetValue (std::size_t idx, std::size_t idy) const
 
G4double GetX (std::size_t index) const
 
G4double GetY (std::size_t index) const
 
G4bool operator!= (const G4Physics2DVector &right) const =delete
 
G4Physics2DVectoroperator= (const G4Physics2DVector &)
 
G4bool operator== (const G4Physics2DVector &right) const =delete
 
void PutValue (std::size_t idx, std::size_t idy, G4double value)
 
void PutVectors (const std::vector< G4double > &vecX, const std::vector< G4double > &vecY)
 
void PutX (std::size_t idx, G4double value)
 
void PutY (std::size_t idy, G4double value)
 
G4bool Retrieve (std::ifstream &fIn)
 
void ScaleVector (G4double factor)
 
void SetBicubicInterpolation (G4bool)
 
void SetVerboseLevel (G4int value)
 
void Store (std::ofstream &fOut) const
 
G4double Value (G4double x, G4double y) const
 
G4double Value (G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
 
 ~G4Physics2DVector ()
 

Protected Member Functions

G4double BicubicInterpolation (const G4double x, const G4double y, const std::size_t idx, const std::size_t idy) const
 
void ClearVectors ()
 
void CopyData (const G4Physics2DVector &vec)
 
std::size_t FindBin (const G4double z, const G4PV2DDataVector &, const std::size_t idz, const std::size_t idzmax) const
 
void PrepareVectors ()
 

Private Member Functions

G4double DerivativeX (std::size_t idx, std::size_t idy, G4double fac) const
 
G4double DerivativeXY (std::size_t idx, std::size_t idy, G4double fac) const
 
G4double DerivativeY (std::size_t idx, std::size_t idy, G4double fac) const
 
G4double InterpolateLinearX (G4PV2DDataVector &v, G4double rand) const
 

Private Attributes

std::size_t numberOfXNodes = 0
 
std::size_t numberOfYNodes = 0
 
G4PhysicsVectorType type = T_G4PhysicsFreeVector
 
G4bool useBicubic = false
 
std::vector< G4PV2DDataVector * > value
 
G4int verboseLevel = 0
 
G4PV2DDataVector xVector
 
G4PV2DDataVector yVector
 

Detailed Description

Definition at line 47 of file G4Physics2DVector.hh.

Constructor & Destructor Documentation

◆ G4Physics2DVector() [1/3]

G4Physics2DVector::G4Physics2DVector ( )

Definition at line 37 of file G4Physics2DVector.cc.

References PrepareVectors().

◆ G4Physics2DVector() [2/3]

G4Physics2DVector::G4Physics2DVector ( std::size_t  nx,
std::size_t  ny 
)
explicit

Definition at line 41 of file G4Physics2DVector.cc.

42{
43 if(nx < 2 || ny < 2)
44 {
46 ed << "G4Physics2DVector is too short: nx= " << nx << " numy= " << ny;
47 G4Exception("G4Physics2DVector::G4Physics2DVector()", "glob03",
48 FatalException, ed, "Both lengths should be above 1");
49 }
50 numberOfXNodes = nx;
51 numberOfYNodes = ny;
53}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::size_t numberOfXNodes
std::size_t numberOfYNodes

References FatalException, G4Exception(), numberOfXNodes, numberOfYNodes, and PrepareVectors().

◆ G4Physics2DVector() [3/3]

G4Physics2DVector::G4Physics2DVector ( const G4Physics2DVector right)

Definition at line 61 of file G4Physics2DVector.cc.

62{
63 type = right.type;
64
67
69 useBicubic = right.useBicubic;
70
71 xVector = right.xVector;
72 yVector = right.yVector;
73
75 CopyData(right);
76}
void CopyData(const G4Physics2DVector &vec)
G4PV2DDataVector xVector
G4PV2DDataVector yVector
G4PhysicsVectorType type

References CopyData(), numberOfXNodes, numberOfYNodes, PrepareVectors(), type, useBicubic, verboseLevel, xVector, and yVector.

◆ ~G4Physics2DVector()

G4Physics2DVector::~G4Physics2DVector ( )

Definition at line 57 of file G4Physics2DVector.cc.

References ClearVectors().

Member Function Documentation

◆ BicubicInterpolation()

G4double G4Physics2DVector::BicubicInterpolation ( const G4double  x,
const G4double  y,
const std::size_t  idx,
const std::size_t  idy 
) const
protected

Definition at line 182 of file G4Physics2DVector.cc.

186{
187 // Bicubic interpolation according to
188 // 1. H.M. Antia, "Numerical Methods for Scientists and Engineers",
189 // MGH, 1991.
190 // 2. W.H. Press et al., "Numerical recipes. The Art of Scientific
191 // Computing", Cambridge University Press, 2007.
192 const G4double x1 = xVector[idx];
193 const G4double x2 = xVector[idx + 1];
194 const G4double y1 = yVector[idy];
195 const G4double y2 = yVector[idy + 1];
196 const G4double f1 = GetValue(idx, idy);
197 const G4double f2 = GetValue(idx + 1, idy);
198 const G4double f3 = GetValue(idx + 1, idy + 1);
199 const G4double f4 = GetValue(idx, idy + 1);
200
201 const G4double dx = x2 - x1;
202 const G4double dy = y2 - y1;
203
204 const G4double h1 = (x - x1) / dx;
205 const G4double h2 = (y - y1) / dy;
206
207 const G4double h12 = h1 * h1;
208 const G4double h13 = h12 * h1;
209 const G4double h22 = h2 * h2;
210 const G4double h23 = h22 * h2;
211
212 // Three derivatives at each of four points (1-4) defining the
213 // subregion are computed by numerical centered differencing from
214 // the functional values already tabulated on the grid.
215
216 const G4double f1x = DerivativeX(idx, idy, dx);
217 const G4double f2x = DerivativeX(idx + 1, idy, dx);
218 const G4double f3x = DerivativeX(idx + 1, idy + 1, dx);
219 const G4double f4x = DerivativeX(idx, idy + 1, dx);
220
221 const G4double f1y = DerivativeY(idx, idy, dy);
222 const G4double f2y = DerivativeY(idx + 1, idy, dy);
223 const G4double f3y = DerivativeY(idx + 1, idy + 1, dy);
224 const G4double f4y = DerivativeY(idx, idy + 1, dy);
225
226 const G4double dxy = dx * dy;
227 const G4double f1xy = DerivativeXY(idx, idy, dxy);
228 const G4double f2xy = DerivativeXY(idx + 1, idy, dxy);
229 const G4double f3xy = DerivativeXY(idx + 1, idy + 1, dxy);
230 const G4double f4xy = DerivativeXY(idx, idy + 1, dxy);
231
232 return f1 + f1y * h2 + (3 * (f4 - f1) - 2 * f1y - f4y) * h22 +
233 (2 * (f1 - f4) + f1y + f4y) * h23 + f1x * h1 + f1xy * h1 * h2 +
234 (3 * (f4x - f1x) - 2 * f1xy - f4xy) * h1 * h22 +
235 (2 * (f1x - f4x) + f1xy + f4xy) * h1 * h23 +
236 (3 * (f2 - f1) - 2 * f1x - f2x) * h12 +
237 (3 * f2y - 3 * f1y - 2 * f1xy - f2xy) * h12 * h2 +
238 (9 * (f1 - f2 + f3 - f4) + 6 * f1x + 3 * f2x - 3 * f3x - 6 * f4x +
239 6 * f1y - 6 * f2y - 3 * f3y + 3 * f4y + 4 * f1xy + 2 * f2xy + f3xy +
240 2 * f4xy) *
241 h12 * h22 +
242 (6 * (-f1 + f2 - f3 + f4) - 4 * f1x - 2 * f2x + 2 * f3x + 4 * f4x -
243 3 * f1y + 3 * f2y + 3 * f3y - 3 * f4y - 2 * f1xy - f2xy - f3xy -
244 2 * f4xy) *
245 h12 * h23 +
246 (2 * (f1 - f2) + f1x + f2x) * h13 +
247 (2 * (f1y - f2y) + f1xy + f2xy) * h13 * h2 +
248 (6 * (-f1 + f2 - f3 + f4) + 3 * (-f1x - f2x + f3x + f4x) - 4 * f1y +
249 4 * f2y + 2 * f3y - 2 * f4y - 2 * f1xy - 2 * f2xy - f3xy - f4xy) *
250 h13 * h22 +
251 (4 * (f1 - f2 + f3 - f4) + 2 * (f1x + f2x - f3x - f4x) +
252 2 * (f1y - f2y - f3y + f4y) + f1xy + f2xy + f3xy + f4xy) *
253 h13 * h23;
254}
double G4double
Definition: G4Types.hh:83
G4double DerivativeX(std::size_t idx, std::size_t idy, G4double fac) const
G4double DerivativeXY(std::size_t idx, std::size_t idy, G4double fac) const
G4double DerivativeY(std::size_t idx, std::size_t idy, G4double fac) const
G4double GetValue(std::size_t idx, std::size_t idy) const

References DerivativeX(), DerivativeXY(), DerivativeY(), GetValue(), xVector, and yVector.

Referenced by Value().

◆ ClearVectors()

void G4Physics2DVector::ClearVectors ( )
protected

Definition at line 117 of file G4Physics2DVector.cc.

118{
119 for(std::size_t j = 0; j < numberOfYNodes; ++j)
120 {
121 delete value[j];
122 }
123}
std::vector< G4PV2DDataVector * > value

References numberOfYNodes, and value.

Referenced by operator=(), PutVectors(), Retrieve(), and ~G4Physics2DVector().

◆ CopyData()

void G4Physics2DVector::CopyData ( const G4Physics2DVector vec)
protected

Definition at line 127 of file G4Physics2DVector.cc.

128{
129 for(std::size_t i = 0; i < numberOfXNodes; ++i)
130 {
131 xVector[i] = right.xVector[i];
132 }
133 for(std::size_t j = 0; j < numberOfYNodes; ++j)
134 {
135 yVector[j] = right.yVector[j];
136 G4PV2DDataVector* v0 = right.value[j];
137 for(std::size_t i = 0; i < numberOfXNodes; ++i)
138 {
139 PutValue(i, j, (*v0)[i]);
140 }
141 }
142}
std::vector< G4double > G4PV2DDataVector
void PutValue(std::size_t idx, std::size_t idy, G4double value)

References numberOfXNodes, numberOfYNodes, PutValue(), value, xVector, and yVector.

Referenced by G4Physics2DVector(), and operator=().

◆ DerivativeX()

G4double G4Physics2DVector::DerivativeX ( std::size_t  idx,
std::size_t  idy,
G4double  fac 
) const
inlineprivate

Referenced by BicubicInterpolation().

◆ DerivativeXY()

G4double G4Physics2DVector::DerivativeXY ( std::size_t  idx,
std::size_t  idy,
G4double  fac 
) const
inlineprivate

Referenced by BicubicInterpolation().

◆ DerivativeY()

G4double G4Physics2DVector::DerivativeY ( std::size_t  idx,
std::size_t  idy,
G4double  fac 
) const
inlineprivate

Referenced by BicubicInterpolation().

◆ FindBin()

std::size_t G4Physics2DVector::FindBin ( const G4double  z,
const G4PV2DDataVector ,
const std::size_t  idz,
const std::size_t  idzmax 
) const
inlineprotected

◆ FindBinLocationX()

std::size_t G4Physics2DVector::FindBinLocationX ( const G4double  x,
const std::size_t  lastidx 
) const
inline

Referenced by Value().

◆ FindBinLocationY()

std::size_t G4Physics2DVector::FindBinLocationY ( const G4double  y,
const std::size_t  lastidy 
) const
inline

Referenced by FindLinearX(), and Value().

◆ FindLinearX() [1/2]

G4double G4Physics2DVector::FindLinearX ( G4double  rand,
G4double  y 
) const
inline

◆ FindLinearX() [2/2]

G4double G4Physics2DVector::FindLinearX ( G4double  rand,
G4double  y,
std::size_t &  lastidy 
) const

Definition at line 388 of file G4Physics2DVector.cc.

390{
392
393 // find bins
394 idy = FindBinLocationY(y, idy);
395
396 G4double x1 = InterpolateLinearX(*(value[idy]), rand);
397 G4double x2 = InterpolateLinearX(*(value[idy + 1]), rand);
398 G4double res = x1;
399 G4double del = yVector[idy + 1] - yVector[idy];
400 if(del != 0.0)
401 {
402 res += (x2 - x1) * (y - yVector[idy]) / del;
403 }
404 return res;
405}
std::size_t FindBinLocationY(const G4double y, const std::size_t lastidy) const
G4double InterpolateLinearX(G4PV2DDataVector &v, G4double rand) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References FindBinLocationY(), InterpolateLinearX(), G4INCL::Math::max(), G4INCL::Math::min(), numberOfYNodes, value, and yVector.

Referenced by G4MuPairProductionModel::FindScaledEnergy().

◆ GetLengthX()

std::size_t G4Physics2DVector::GetLengthX ( ) const
inline

◆ GetLengthY()

std::size_t G4Physics2DVector::GetLengthY ( ) const
inline

◆ GetType()

G4PhysicsVectorType G4Physics2DVector::GetType ( ) const
inline

◆ GetValue()

G4double G4Physics2DVector::GetValue ( std::size_t  idx,
std::size_t  idy 
) const
inline

◆ GetX()

G4double G4Physics2DVector::GetX ( std::size_t  index) const
inline

◆ GetY()

G4double G4Physics2DVector::GetY ( std::size_t  index) const
inline

◆ InterpolateLinearX()

G4double G4Physics2DVector::InterpolateLinearX ( G4PV2DDataVector v,
G4double  rand 
) const
private

Definition at line 409 of file G4Physics2DVector.cc.

411{
412 std::size_t nn = v.size();
413 if(1 >= nn)
414 {
415 return 0.0;
416 }
417 std::size_t n1 = 0;
418 std::size_t n2 = nn / 2;
419 std::size_t n3 = nn - 1;
420 G4double y = rand * v[n3];
421 while(n1 + 1 != n3)
422 {
423 if(y > v[n2])
424 {
425 n1 = n2;
426 }
427 else
428 {
429 n3 = n2;
430 }
431 n2 = (n3 + n1 + 1) / 2;
432 }
433 G4double res = xVector[n1];
434 G4double del = v[n3] - v[n1];
435 if(del > 0.0)
436 {
437 res += (y - v[n1]) * (xVector[n3] - res) / del;
438 }
439 return res;
440}

References G4InuclParticleNames::nn, and xVector.

Referenced by FindLinearX().

◆ operator!=()

G4bool G4Physics2DVector::operator!= ( const G4Physics2DVector right) const
delete

◆ operator=()

G4Physics2DVector & G4Physics2DVector::operator= ( const G4Physics2DVector right)

Definition at line 80 of file G4Physics2DVector.cc.

81{
82 if(&right == this)
83 {
84 return *this;
85 }
87
88 type = right.type;
89
92
94 useBicubic = right.useBicubic;
95
97 CopyData(right);
98
99 return *this;
100}

References ClearVectors(), CopyData(), numberOfXNodes, numberOfYNodes, PrepareVectors(), type, useBicubic, and verboseLevel.

◆ operator==()

G4bool G4Physics2DVector::operator== ( const G4Physics2DVector right) const
delete

◆ PrepareVectors()

void G4Physics2DVector::PrepareVectors ( )
protected

Definition at line 104 of file G4Physics2DVector.cc.

105{
106 xVector.resize(numberOfXNodes, 0.);
107 yVector.resize(numberOfYNodes, 0.);
108 value.resize(numberOfYNodes);
109 for(std::size_t j = 0; j < numberOfYNodes; ++j)
110 {
112 }
113}

References numberOfXNodes, numberOfYNodes, value, xVector, and yVector.

Referenced by G4Physics2DVector(), operator=(), PutVectors(), and Retrieve().

◆ PutValue()

void G4Physics2DVector::PutValue ( std::size_t  idx,
std::size_t  idy,
G4double  value 
)
inline

◆ PutVectors()

void G4Physics2DVector::PutVectors ( const std::vector< G4double > &  vecX,
const std::vector< G4double > &  vecY 
)

Definition at line 258 of file G4Physics2DVector.cc.

260{
261 ClearVectors();
262 std::size_t nx = vecX.size();
263 std::size_t ny = vecY.size();
264 if(nx < 2 || ny < 2)
265 {
267 ed << "G4Physics2DVector is too short: nx= " << nx << " ny= " << ny;
268 G4Exception("G4Physics2DVector::PutVectors()", "glob03", FatalException, ed,
269 "Both lengths should be above 1");
270 }
271
272 numberOfXNodes = nx;
273 numberOfYNodes = ny;
275 for(std::size_t i = 0; i < nx; ++i)
276 {
277 xVector[i] = vecX[i];
278 }
279 for(std::size_t j = 0; j < ny; ++j)
280 {
281 yVector[j] = vecY[j];
282 }
283}

References ClearVectors(), FatalException, G4Exception(), numberOfXNodes, numberOfYNodes, PrepareVectors(), xVector, and yVector.

◆ PutX()

void G4Physics2DVector::PutX ( std::size_t  idx,
G4double  value 
)
inline

◆ PutY()

void G4Physics2DVector::PutY ( std::size_t  idy,
G4double  value 
)
inline

◆ Retrieve()

G4bool G4Physics2DVector::Retrieve ( std::ifstream &  fIn)

Definition at line 320 of file G4Physics2DVector.cc.

321{
322 // initialisation
323 ClearVectors();
324
325 // binning
326 G4int k, nx, ny;
327 in >> k >> nx >> ny;
328 if(in.fail() || 2 > nx || 2 > ny || nx >= INT_MAX || ny >= INT_MAX)
329 {
330 return false;
331 }
332 numberOfXNodes = nx;
333 numberOfYNodes = ny;
336
337 // contents
338 G4double val;
339 for(G4int i = 0; i < nx; ++i)
340 {
341 in >> xVector[i];
342 if(in.fail())
343 {
344 return false;
345 }
346 }
347 for(G4int j = 0; j < ny; ++j)
348 {
349 in >> yVector[j];
350 if(in.fail())
351 {
352 return false;
353 }
354 }
355 for(G4int j = 0; j < ny; ++j)
356 {
357 for(G4int i = 0; i < nx; ++i)
358 {
359 in >> val;
360 if(in.fail())
361 {
362 return false;
363 }
364 PutValue(i, j, val);
365 }
366 }
367 in.close();
368 return true;
369}
G4PhysicsVectorType
int G4int
Definition: G4Types.hh:85
#define INT_MAX
Definition: templates.hh:90

References ClearVectors(), INT_MAX, numberOfXNodes, numberOfYNodes, PrepareVectors(), PutValue(), type, xVector, and yVector.

Referenced by G4SeltzerBergerModel::ReadData(), G4LivermoreBremsstrahlungModel::ReadData(), G4OpticalSurface::ReadDichroicFile(), and G4MuPairProductionModel::RetrieveTables().

◆ ScaleVector()

void G4Physics2DVector::ScaleVector ( G4double  factor)

Definition at line 373 of file G4Physics2DVector.cc.

374{
375 G4double val;
376 for(std::size_t j = 0; j < numberOfYNodes; ++j)
377 {
378 for(std::size_t i = 0; i < numberOfXNodes; ++i)
379 {
380 val = GetValue(i, j) * factor;
381 PutValue(i, j, val);
382 }
383 }
384}

References GetValue(), numberOfXNodes, numberOfYNodes, and PutValue().

◆ SetBicubicInterpolation()

void G4Physics2DVector::SetBicubicInterpolation ( G4bool  )
inline

◆ SetVerboseLevel()

void G4Physics2DVector::SetVerboseLevel ( G4int  value)
inline

◆ Store()

void G4Physics2DVector::Store ( std::ofstream &  fOut) const

Definition at line 287 of file G4Physics2DVector.cc.

288{
289 // binning
290 G4int prec = out.precision();
291 out << G4int(type) << " " << numberOfXNodes << " " << numberOfYNodes
292 << G4endl;
293 out << std::setprecision(8);
294
295 // contents
296 for(std::size_t i = 0; i < numberOfXNodes - 1; ++i)
297 {
298 out << xVector[i] << " ";
299 }
300 out << xVector[numberOfXNodes - 1] << G4endl;
301 for(std::size_t j = 0; j < numberOfYNodes - 1; ++j)
302 {
303 out << yVector[j] << " ";
304 }
305 out << yVector[numberOfYNodes - 1] << G4endl;
306 for(std::size_t j = 0; j < numberOfYNodes; ++j)
307 {
308 for(std::size_t i = 0; i < numberOfXNodes - 1; ++i)
309 {
310 out << GetValue(i, j) << " ";
311 }
312 out << GetValue(numberOfXNodes - 1, j) << G4endl;
313 }
314 out.precision(prec);
315 out.close();
316}
#define G4endl
Definition: G4ios.hh:57
static const double prec
Definition: RanecuEngine.cc:61

References G4endl, GetValue(), numberOfXNodes, numberOfYNodes, CLHEP::prec, type, xVector, and yVector.

Referenced by G4MuPairProductionModel::StoreTables().

◆ Value() [1/2]

G4double G4Physics2DVector::Value ( G4double  x,
G4double  y 
) const

◆ Value() [2/2]

G4double G4Physics2DVector::Value ( G4double  x,
G4double  y,
std::size_t &  lastidx,
std::size_t &  lastidy 
) const

Definition at line 146 of file G4Physics2DVector.cc.

148{
149 // no interpolation outside the table
150 const G4double x =
152 const G4double y =
154
155 // find bins
156 idx = FindBinLocationX(x, idx);
157 idy = FindBinLocationY(y, idy);
158
159 // interpolate
160 if(useBicubic)
161 {
162 return BicubicInterpolation(x, y, idx, idy);
163 }
164 else
165 {
166 const G4double x1 = xVector[idx];
167 const G4double x2 = xVector[idx + 1];
168 const G4double y1 = yVector[idy];
169 const G4double y2 = yVector[idy + 1];
170 const G4double v11 = GetValue(idx, idy);
171 const G4double v12 = GetValue(idx + 1, idy);
172 const G4double v21 = GetValue(idx, idy + 1);
173 const G4double v22 = GetValue(idx + 1, idy + 1);
174 return ((y2 - y) * (v11 * (x2 - x) + v12 * (x - x1)) +
175 ((y - y1) * (v21 * (x2 - x) + v22 * (x - x1)))) /
176 ((x2 - x1) * (y2 - y1));
177 }
178}
std::size_t FindBinLocationX(const G4double x, const std::size_t lastidx) const
G4double BicubicInterpolation(const G4double x, const G4double y, const std::size_t idx, const std::size_t idy) const

References BicubicInterpolation(), FindBinLocationX(), FindBinLocationY(), GetValue(), G4INCL::Math::max(), G4INCL::Math::min(), numberOfXNodes, numberOfYNodes, useBicubic, xVector, and yVector.

Referenced by G4eDPWAElasticDCS::ComputeCSPerAtom(), G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom(), G4SeltzerBergerModel::ComputeDXSectionPerAtom(), G4OpBoundaryProcess::DielectricDichroic(), G4MuPairProductionModel::FindScaledEnergy(), G4ChannelingECHARM::GetEC(), G4eDPWAElasticDCS::LoadDCSForZ(), G4SeltzerBergerModel::ReadData(), G4LivermoreBremsstrahlungModel::ReadData(), G4SeltzerBergerModel::SampleEnergyTransfer(), and G4LivermoreBremsstrahlungModel::SampleSecondaries().

Field Documentation

◆ numberOfXNodes

std::size_t G4Physics2DVector::numberOfXNodes = 0
private

◆ numberOfYNodes

std::size_t G4Physics2DVector::numberOfYNodes = 0
private

◆ type

G4PhysicsVectorType G4Physics2DVector::type = T_G4PhysicsFreeVector
private

Definition at line 145 of file G4Physics2DVector.hh.

Referenced by G4Physics2DVector(), operator=(), Retrieve(), and Store().

◆ useBicubic

G4bool G4Physics2DVector::useBicubic = false
private

Definition at line 156 of file G4Physics2DVector.hh.

Referenced by G4Physics2DVector(), operator=(), and Value().

◆ value

std::vector<G4PV2DDataVector*> G4Physics2DVector::value
private

Definition at line 153 of file G4Physics2DVector.hh.

Referenced by ClearVectors(), CopyData(), FindLinearX(), and PrepareVectors().

◆ verboseLevel

G4int G4Physics2DVector::verboseLevel = 0
private

Definition at line 155 of file G4Physics2DVector.hh.

Referenced by G4Physics2DVector(), and operator=().

◆ xVector

G4PV2DDataVector G4Physics2DVector::xVector
private

◆ yVector

G4PV2DDataVector G4Physics2DVector::yVector
private

The documentation for this class was generated from the following files: