QwGeant4
QweakSimFieldMap< value_t, value_n > Class Template Reference

A multi-dimensional grid of values with interpolation methods. More...

#include <QweakSimFieldMap.hh>

Public Types

enum  debug_t { kError, kWarning, kVerbose, kDebug }
 Debug level. More...
 

Public Member Functions

 QweakSimFieldMap (const unsigned int ndim=1)
 Constructor with number of dimensions. More...
 
 QweakSimFieldMap (const std::vector< coord_t > &min, const std::vector< coord_t > &max, const std::vector< coord_t > &step)
 Constructor with minimum, maximum, and step size. More...
 
 QweakSimFieldMap (const std::string &filename)
 Constructor with file name. More...
 
virtual ~QweakSimFieldMap ()
 Destructor. More...
 
void SetDebugLevel (debug_t debug)
 Set debug level. More...
 
void SetDimensions (const unsigned int ndim)
 Set the number of coordinate dimensions and resize vectors. More...
 
void SetMinimumMaximumStep (const coord_t min, const coord_t max, const coord_t step)
 Set minimum, maximum, and step size to single values. More...
 
void SetMinimumMaximumStep (const std::vector< coord_t > &min, const std::vector< coord_t > &max, const std::vector< coord_t > &step)
 Set minimum, maximum, and step size to different values. More...
 
coord_t GetMinimum (const unsigned int dim) const
 Get minimum in dimension. More...
 
coord_t GetMaximum (const unsigned int dim) const
 Get maximum in dimension. More...
 
coord_t GetStepSize (const unsigned int dim) const
 Get minimum in dimension. More...
 
unsigned int GetMaximumEntries () const
 Get the maximum number of entries. More...
 
unsigned int GetCurrentEntries () const
 Get the current number of entries. More...
 
unsigned int GetWrapCoordinate (const unsigned int dim) const
 Get wrapping coordinate. More...
 
void SetWrapCoordinate (const unsigned int dim, const size_t wrap=1)
 Set wrapping coordinate. More...
 
void SetWrapCoordinate (const std::vector< size_t > &wrap)
 
int GetDataReductionFactor (const unsigned int dim) const
 Get data reduction factor. More...
 
void SetDataReductionFactor (const unsigned int dim, const unsigned int redux)
 Set data reduction factor. More...
 
void SetDataReductionFactor (const unsigned int redux)
 
void SetDataReductionFactor (const std::vector< unsigned int > &redux)
 
void SetInterpolationMethod (const EInterpolationMethod method)
 Set the interpolation method. More...
 
EInterpolationMethod GetInterpolationMethod () const
 Get the interpolation method. More...
 
void PrintCoverage (const unsigned int dim)
 Print coverage map for all bins in one dimension. More...
 
bool InBounds (const coord_t *coord) const
 Return true if the coordinate is in bounds. More...
 
Functions to write grid values
bool Set (const coord_t &coord, const value_t &value)
 Set a single value at a coordinate (false if not possible) More...
 
bool Set (const coord_t *coord, const value_t &value)
 Set a single value at a coordinate (false if not possible) More...
 
bool Set (const coord_t *coord, const value_t *value)
 Set a set of values at a coordinate (false if not possible) More...
 
bool Set (const unsigned int linear_index, const value_t &value)
 Set a single value at a linearized index (false if not possible) More...
 
bool Set (const unsigned int linear_index, const value_t *value)
 Set a set of values at a linearized index (false if not possible) More...
 
bool Set (const unsigned int *cell_index, const value_t &value)
 Set a single value at a grid point (false if out of bounds) More...
 
bool Set (const unsigned int *cell_index, const value_t *value)
 Set a set of values at a grid point (false if out of bounds) More...
 
Functions to retrieve interpolated values
value_t GetValue (const coord_t &coord) const
 Get the interpolated value at coordinate (zero if out of bounds) More...
 
value_t GetValue (const coord_t *coord) const
 Get the interpolated value at coordinate (zero if out of bounds) More...
 
bool GetValue (const coord_t *coord, double &value) const
 Get the interpolated value at coordinate (zero if out of bounds) More...
 
bool GetValue (const coord_t *coord, double *value) const
 Get the interpolated value at coordinate (zero if out of bounds) More...
 
File reading and writing
bool WriteText (std::ostream &stream) const
 Write the grid as text. More...
 
bool WriteTextFile (const std::string &filename) const
 Write the grid to text file. More...
 
bool WriteTextScreen () const
 Write the grid to screen. More...
 
bool ReadText (std::istream &stream)
 Read the grid from text. More...
 
bool ReadTextFile (const std::string &filename)
 Read the grid from text file. More...
 
bool WriteBinaryFile (const std::string &filename) const
 Write the grid values to binary file. More...
 
bool ReadBinaryFile (const std::string &filename)
 Read the grid values from binary file. More...
 
Indexing functions (publicly available and unchecked)
unsigned int Index (const coord_t *coord) const
 Return the linearized index based on the point coordinates (unchecked) More...
 
unsigned int Index (const unsigned int *cell_index) const
 Return the linearized index based on the cell indices (unchecked) More...
 
unsigned int Index (const unsigned int *cell_index, const unsigned int offset) const
 Return the linearized index based on the cell indices and offset (unchecked) More...
 
void Cell (const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
 Return the cell index and local coordinates in one dimension (unchecked) More...
 
void Cell (const coord_t *coord, unsigned int *cell_index, double *cell_local) const
 Return the cell index and local coordinates (unchecked) More...
 
void Cell (const coord_t *coord, unsigned int *cell_index) const
 Return the cell index (unchecked) More...
 
void Cell (unsigned int linear_index, unsigned int *cell_index) const
 Return the cell index based on the linearized index (unchecked) More...
 
void Coord (const unsigned int *cell_index, coord_t *coord) const
 Return the coordinates based on the cell index (unchecked) More...
 
void Coord (const unsigned int linear_index, coord_t *coord) const
 Return the coordinates based on the linearized index (unchecked) More...
 

Data Fields

enum QweakSimFieldMap::debug_t fDebug
 

Private Member Functions

void Nearest (const coord_t *coord, unsigned int *cell_index) const
 Return the cell index closest to the coordinate (could be above) (unchecked) More...
 
bool Linear (const coord_t *coord, value_t *value) const
 Linear interpolation (unchecked) More...
 
bool NearestNeighbor (const coord_t *coord, value_t *value) const
 Nearest-neighbor (unchecked) More...
 
bool Check (const coord_t *coord) const
 Check for boundaries with coordinate argument. More...
 
bool Check (const unsigned int *cell_index) const
 Check for boundaries with cell index argument. More...
 
bool Check (const unsigned int linear_index) const
 Check for boundaries with linearized index argument. More...
 
value_t Get (const unsigned int *cell_index) const
 Get a single value by cell index (unchecked) More...
 
value_t Get (const unsigned int index) const
 Get a single value by linearized index (unchecked) More...
 
bool Get (const unsigned int index, value_t *value) const
 Get a vector value by linearized index (unchecked) More...
 

Private Attributes

unsigned int fNDim
 Number of dimensions in coordinates. More...
 
std::vector< coord_tfMin
 Minimum in each dimension. More...
 
std::vector< coord_tfMax
 Maximum in each dimension. More...
 
std::vector< coord_tfStep
 Step size in each dimension. More...
 
std::vector< size_t > fSize
 Number of points in each dimension. More...
 
std::vector< size_t > fWrap
 Wrap around this coordinate. More...
 
std::vector< size_t > fRedux
 Data reduction factor. More...
 
std::vector< size_t > fExtent
 
std::vector< value_t > fValues [value_n]
 Table with pointers to arrays of values. More...
 
unsigned int fCurrentEntries
 Number of values read in. More...
 
unsigned int fMaximumEntries
 Maximum number of values. More...
 
EInterpolationMethod fInterpolationMethod
 Interpolation method. More...
 

Detailed Description

template<class value_t = float, unsigned int value_n = 1>
class QweakSimFieldMap< value_t, value_n >

A multi-dimensional grid of values with interpolation methods.

This class provides various interpolation methods on a multi-dimensional grid of multi-dimensional values. Linear interpolation and nearest-neighbor are implemented for all dimensions.

The template arguments are the internal storage data type (defaults to float) and the number of dimensions of the stored data at each grid point (defaults to one). The dimension of the grid itself is set through the constructor. To describe a double vector field with 5 components on a 3-dimensional grid, you would write

* QwInterpolate<float,5> *field = new QweakSimFieldMap<float,5>(3);
*

The minimum, maximum, and step size of the grid have to be known before the values are filled.

Definition at line 51 of file QweakSimFieldMap.hh.

Member Enumeration Documentation

template<class value_t = float, unsigned int value_n = 1>
enum QweakSimFieldMap::debug_t

Debug level.

Enumerator
kError 
kWarning 
kVerbose 
kDebug 

Definition at line 78 of file QweakSimFieldMap.hh.

Constructor & Destructor Documentation

template<class value_t = float, unsigned int value_n = 1>
QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap ( const unsigned int  ndim = 1)
inline

Constructor with number of dimensions.

Definition at line 56 of file QweakSimFieldMap.hh.

References kMultiLinear, QweakSimFieldMap< value_t, value_n >::SetDimensions(), and QweakSimFieldMap< value_t, value_n >::SetInterpolationMethod().

56  {
57  SetDimensions(ndim);
59  };
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
void SetInterpolationMethod(const EInterpolationMethod method)
Set the interpolation method.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap ( const std::vector< coord_t > &  min,
const std::vector< coord_t > &  max,
const std::vector< coord_t > &  step 
)
inline

Constructor with minimum, maximum, and step size.

Definition at line 61 of file QweakSimFieldMap.hh.

References kMultiLinear, QweakSimFieldMap< value_t, value_n >::SetDimensions(), QweakSimFieldMap< value_t, value_n >::SetInterpolationMethod(), and QweakSimFieldMap< value_t, value_n >::SetMinimumMaximumStep().

63  {
64  SetDimensions(min.size());
66  SetMinimumMaximumStep(min,max,step);
67  };
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
void SetInterpolationMethod(const EInterpolationMethod method)
Set the interpolation method.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap ( const std::string &  filename)
inline

Constructor with file name.

Definition at line 69 of file QweakSimFieldMap.hh.

References kMultiLinear, QweakSimFieldMap< value_t, value_n >::ReadBinaryFile(), and QweakSimFieldMap< value_t, value_n >::SetInterpolationMethod().

69  {
70  ReadBinaryFile(filename);
72  };
bool ReadBinaryFile(const std::string &filename)
Read the grid values from binary file.
void SetInterpolationMethod(const EInterpolationMethod method)
Set the interpolation method.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
virtual QweakSimFieldMap< value_t, value_n >::~QweakSimFieldMap ( )
inlinevirtual

Destructor.

Definition at line 74 of file QweakSimFieldMap.hh.

74 { };

Member Function Documentation

template<class value_t , unsigned int value_n>
void QweakSimFieldMap< value_t, value_n >::Cell ( const coord_t  coord,
unsigned int &  cell_index,
double &  cell_local,
const unsigned int  dim 
) const
inline

Return the cell index and local coordinates in one dimension (unchecked)

Return the cell index and local coordinates in one dimension (unchecked)

Parameters
coordPoint coordinate in one dimension
cell_indexCell index of the point (reference)
cell_localLocal coordinates in cell (reference)
dimDimension

Definition at line 723 of file QweakSimFieldMap.hh.

Referenced by QweakSimFieldMap< value_t, value_n >::Index(), QweakSimFieldMap< value_t, value_n >::Nearest(), and QweakSimFieldMap< value_t, value_n >::PrintCoverage().

728 {
729  // Normalize coordinate (positive, with integers on grid points)
730  double norm_coord = (coord - fMin[dim]) / fStep[dim];
731  // Split in fractional and integer part
732  double frac_part, int_part;
733  frac_part = modf(norm_coord, &int_part);
734  cell_local = frac_part;
735  cell_index = static_cast<int>(int_part); // cast to integer
736  // Wrap index
737  if (fWrap[dim] > 0) cell_index %= (fSize[dim] - fWrap[dim]);
738 }
std::vector< coord_t > fStep
Step size in each dimension.
std::vector< size_t > fWrap
Wrap around this coordinate.
std::vector< size_t > fSize
Number of points in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
void QweakSimFieldMap< value_t, value_n >::Cell ( const coord_t coord,
unsigned int *  cell_index,
double *  cell_local 
) const
inline

Return the cell index and local coordinates (unchecked)

Return the cell index and local coordinates (unchecked)

Parameters
coordPoint coordinates
cell_indexCell index of the point (reference)
cell_localLocal coordinates in cell (reference)

Definition at line 747 of file QweakSimFieldMap.hh.

751 {
752  // Get cell index and local coordinates in each dimension
753  for (unsigned int dim = 0; dim < fNDim; dim++)
754  Cell(coord[dim], cell_index[dim], cell_local[dim], dim);
755 }
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
template<class value_t , unsigned int value_n>
void QweakSimFieldMap< value_t, value_n >::Cell ( const coord_t coord,
unsigned int *  cell_index 
) const
inline

Return the cell index (unchecked)

Return the cell index (unchecked)

Parameters
coordPoint coordinates
cell_indexCell index of the point (reference)

Definition at line 763 of file QweakSimFieldMap.hh.

766 {
767  // Get cell index and ignore local coordinates
768  double* cell_local = new double[fNDim];
769  Cell(coord, cell_index, cell_local);
770  delete[] cell_local;
771 }
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
template<class value_t , unsigned int value_n>
void QweakSimFieldMap< value_t, value_n >::Cell ( unsigned int  linear_index,
unsigned int *  cell_index 
) const
inline

Return the cell index based on the linearized index (unchecked)

Return the cell index based on the linearized index (unchecked)

Parameters
linear_indexLinearized index
cell_indexCell index (reference)

Definition at line 779 of file QweakSimFieldMap.hh.

782 {
783  // This does not work with unsigned int, because that is always >= 0 and wraps
784  for (int dim = fNDim-1; dim >= 0; dim--) {
785  cell_index[dim] = linear_index / fExtent[dim];
786  linear_index -= cell_index[dim] * fExtent[dim];
787  }
788 }
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::Check ( const coord_t coord) const
inlineprivate

Check for boundaries with coordinate argument.

Check whether the point is inside the valid region

Parameters
coordPoint coordinates
Returns
True if the coordinates are in the valid region

Definition at line 619 of file QweakSimFieldMap.hh.

References mycout, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::GetValue(), QweakSimFieldMap< value_t, value_n >::InBounds(), and QweakSimFieldMap< value_t, value_n >::Set().

620 {
621  for (unsigned int dim = 0; dim < fNDim; dim++) {
622  if (fWrap[dim] == 0 && (coord[dim] < fMin[dim] || coord[dim] > fMax[dim])) {
623  if (fDebug >= kDebug) {
624  mycout << "QweakSimFieldMap::Check(coord_t*): "
625  << "coord[" << dim << "] out of bounds "
626  << fMin[dim] << "-" << fMax[dim] << "." << myendl;
627  }
628  return false;
629  }
630  }
631  // Otherwise
632  return true;
633 }
std::vector< coord_t > fMax
Maximum in each dimension.
enum QweakSimFieldMap::debug_t fDebug
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fWrap
Wrap around this coordinate.
std::vector< coord_t > fMin
Minimum in each dimension.

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::Check ( const unsigned int *  cell_index) const
inlineprivate

Check for boundaries with cell index argument.

Check whether the cell index is inside the valid region

Parameters
cell_indexCell index
Returns
True if the cell index is in the valid region

Definition at line 641 of file QweakSimFieldMap.hh.

References mycout, and myendl.

642 {
643  for (unsigned int dim = 0; dim < fNDim; dim++) {
644  if (cell_index[dim] >= fSize[dim]) {
645  if (fDebug >= kDebug) {
646  mycout << "QweakSimFieldMap::Check(unsigned int*): "
647  << "cell index[" << dim << "] " << cell_index[dim] << " larger "
648  << "than maximum " << fSize[dim] << "." << myendl;
649  }
650  return false;
651  }
652  }
653  // Otherwise
654  return true;
655 }
enum QweakSimFieldMap::debug_t fDebug
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fSize
Number of points in each dimension.
template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::Check ( const unsigned int  linear_index) const
inlineprivate

Check for boundaries with linearized index argument.

Check whether the linearized index is inside the valid region

Parameters
linear_indexLinearized index
Returns
True if the cell index is in the valid region

Definition at line 663 of file QweakSimFieldMap.hh.

References mycout, and myendl.

664 {
665  if (linear_index >= fMaximumEntries) {
666  if (fDebug >= kDebug) {
667  mycout << "QweakSimFieldMap::Check(unsigned int): "
668  << "linear index " << linear_index << " larger "
669  << "than maximum " << fMaximumEntries << "." << myendl;
670  }
671  return false;
672  }
673  // Otherwise
674  return true;
675 }
enum QweakSimFieldMap::debug_t fDebug
unsigned int fMaximumEntries
Maximum number of values.
#define myendl
#define mycout
template<class value_t , unsigned int value_n>
void QweakSimFieldMap< value_t, value_n >::Coord ( const unsigned int *  cell_index,
coord_t coord 
) const
inline

Return the coordinates based on the cell index (unchecked)

Return the coordinates based on the cell index (unchecked)

Parameters
cell_indexCell index
coordPoint coordinates (reference)

Definition at line 796 of file QweakSimFieldMap.hh.

Referenced by QweakSimFieldMap< value_t, value_n >::PrintCoverage().

799 {
800  for (unsigned int dim = 0; dim < fNDim; dim++)
801  coord[dim] = fMin[dim] + cell_index[dim] * fStep[dim];
802 }
std::vector< coord_t > fStep
Step size in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< coord_t > fMin
Minimum in each dimension.

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
void QweakSimFieldMap< value_t, value_n >::Coord ( const unsigned int  linear_index,
coord_t coord 
) const
inline

Return the coordinates based on the linearized index (unchecked)

Return the coordinates based on the linearized index (unchecked)

Parameters
linear_indexLinearized index
coordPoint coordinates (reference)

Definition at line 810 of file QweakSimFieldMap.hh.

813 {
814  unsigned int* cell_index = new unsigned int[fNDim];
815  Cell(linear_index,cell_index);
816  Coord(cell_index,coord);
817  delete[] cell_index;
818 }
void Coord(const unsigned int *cell_index, coord_t *coord) const
Return the coordinates based on the cell index (unchecked)
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
template<class value_t = float, unsigned int value_n = 1>
value_t QweakSimFieldMap< value_t, value_n >::Get ( const unsigned int *  cell_index) const
inlineprivate

Get a single value by cell index (unchecked)

Definition at line 526 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fValues, QweakSimFieldMap< value_t, value_n >::Index(), mycerr, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::PrintCoverage().

526  {
527  if (value_n != 1) {
528  mycerr << "QweakSimFieldMap::Get(unsigned int*): "
529  << "only valid for one-dimensional fields." << myendl;
530  return 0;
531  }
532  return fValues[0][Index(cell_index)];
533  };
#define mycerr
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
#define myendl
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
value_t QweakSimFieldMap< value_t, value_n >::Get ( const unsigned int  index) const
inlineprivate

Get a single value by linearized index (unchecked)

Definition at line 536 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fValues, mycerr, and myendl.

536  {
537  if (value_n != 1) {
538  mycerr << "QweakSimFieldMap::Get(unsigned int): "
539  << "only valid for one-dimensional fields." << myendl;
540  return 0;
541  }
542  return fValues[0][index];
543  };
#define mycerr
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
#define myendl
template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Get ( const unsigned int  index,
value_t *  value 
) const
inlineprivate

Get a vector value by linearized index (unchecked)

Definition at line 545 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fValues.

545  {
546  for (unsigned int i = 0; i < value_n; i++)
547  value[i] = fValues[i][index];
548  return true;
549  };
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
template<class value_t = float, unsigned int value_n = 1>
unsigned int QweakSimFieldMap< value_t, value_n >::GetCurrentEntries ( ) const
inline

Get the current number of entries.

Definition at line 172 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fCurrentEntries.

172 { return fCurrentEntries; };
unsigned int fCurrentEntries
Number of values read in.
template<class value_t = float, unsigned int value_n = 1>
int QweakSimFieldMap< value_t, value_n >::GetDataReductionFactor ( const unsigned int  dim) const
inline

Get data reduction factor.

Definition at line 186 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fRedux.

187  { return fRedux.at(dim); }
std::vector< size_t > fRedux
Data reduction factor.
template<class value_t = float, unsigned int value_n = 1>
EInterpolationMethod QweakSimFieldMap< value_t, value_n >::GetInterpolationMethod ( ) const
inline

Get the interpolation method.

Definition at line 204 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fInterpolationMethod.

205  { return fInterpolationMethod; };
EInterpolationMethod fInterpolationMethod
Interpolation method.
template<class value_t = float, unsigned int value_n = 1>
coord_t QweakSimFieldMap< value_t, value_n >::GetMaximum ( const unsigned int  dim) const
inline

Get maximum in dimension.

Definition at line 166 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fMax.

166 { return fMax.at(dim); };
std::vector< coord_t > fMax
Maximum in each dimension.
template<class value_t = float, unsigned int value_n = 1>
unsigned int QweakSimFieldMap< value_t, value_n >::GetMaximumEntries ( ) const
inline

Get the maximum number of entries.

Definition at line 170 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fMaximumEntries.

170 { return fMaximumEntries; };
unsigned int fMaximumEntries
Maximum number of values.
template<class value_t = float, unsigned int value_n = 1>
coord_t QweakSimFieldMap< value_t, value_n >::GetMinimum ( const unsigned int  dim) const
inline

Get minimum in dimension.

Definition at line 164 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fMin.

164 { return fMin.at(dim); };
std::vector< coord_t > fMin
Minimum in each dimension.
template<class value_t = float, unsigned int value_n = 1>
coord_t QweakSimFieldMap< value_t, value_n >::GetStepSize ( const unsigned int  dim) const
inline

Get minimum in dimension.

Definition at line 168 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fStep.

168 { return fStep.at(dim); };
std::vector< coord_t > fStep
Step size in each dimension.
template<class value_t = float, unsigned int value_n = 1>
value_t QweakSimFieldMap< value_t, value_n >::GetValue ( const coord_t coord) const
inline

Get the interpolated value at coordinate (zero if out of bounds)

Definition at line 347 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fNDim, mycerr, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::GetValue().

347  {
348  if (value_n != 1) {
349  mycerr << "QweakSimFieldMap::GetValue(coord_t&): "
350  << "only valid for one-dimensional fields." << myendl;
351  return false;
352  }
353  if (fNDim != 1) {
354  mycerr << "QweakSimFieldMap::GetValue(coord_t&): "
355  << "only valid for one-dimensional grids." << myendl;
356  return false;
357  }
358  value_t value;
359  if (GetValue(&coord, &value)) return value;
360  else return 0; // interpolation failed
361  };
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)
#define mycerr
#define myendl
unsigned int fNDim
Number of dimensions in coordinates.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
value_t QweakSimFieldMap< value_t, value_n >::GetValue ( const coord_t coord) const
inline

Get the interpolated value at coordinate (zero if out of bounds)

Definition at line 363 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::GetValue(), mycerr, and myendl.

363  {
364  if (value_n != 1) {
365  mycerr << "QweakSimFieldMap::GetValue(coord_t*): "
366  << "only valid for one-dimensional fields." << myendl;
367  return false;
368  }
369  value_t value;
370  if (GetValue(coord, &value)) return value;
371  else return 0; // interpolation failed
372  };
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)
#define mycerr
#define myendl

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::GetValue ( const coord_t coord,
double &  value 
) const
inline

Get the interpolated value at coordinate (zero if out of bounds)

Definition at line 374 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::GetValue(), mycerr, and myendl.

374  {
375  if (value_n != 1) {
376  mycerr << "QweakSimFieldMap::GetValue(coord_t*,value_t&): "
377  << "only valid for one-dimensional fields." << myendl;
378  return false;
379  }
380  return GetValue(coord, &value);
381  };
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)
#define mycerr
#define myendl

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::GetValue ( const coord_t coord,
double *  value 
) const
inline

Get the interpolated value at coordinate (zero if out of bounds)

Definition at line 383 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Check(), QweakSimFieldMap< value_t, value_n >::fDebug, QweakSimFieldMap< value_t, value_n >::fInterpolationMethod, QweakSimFieldMap< value_t, value_n >::kDebug, kMultiLinear, kNearestNeighbor, QweakSimFieldMap< value_t, value_n >::Linear(), mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::NearestNeighbor().

383  {
384  for (unsigned int i = 0; i < value_n; i++) value[i] = 0.0; // zero
385  if (! Check(coord)) {
386  if (fDebug >= kDebug) {
387  mycerr << "QweakSimFieldMap::GetValue(coord_t*,value_t*): "
388  << "coordinate is out of boundary." << myendl;
389  }
390  return false; // out of bounds
391  }
392  value_t result[value_n]; // we need a local copy of type value_t
393  switch (fInterpolationMethod) {
394  // return interpolation status
395  case kMultiLinear:
396  if (! Linear(coord, result)) {
397  if (fDebug >= kDebug) {
398  mycerr << "QweakSimFieldMap::GetValue(coord_t*,value_t*): "
399  << "linear interpolation failed." << myendl;
400  }
401  return false;
402  }
403  break;
404  case kNearestNeighbor:
405  if (! NearestNeighbor(coord, result)) {
406  if (fDebug >= kDebug) {
407  mycerr << "QweakSimFieldMap::GetValue(coord_t*,value_t*): "
408  << "nearest-neighbor interpolation failed." << myendl;
409  }
410  return false;
411  }
412  break;
413  default:
414  mycerr << "QweakSimFieldMap::GetValue(coord_t*,value_t*): "
415  << "unknown interpolation method." << myendl;
416  return false;
417  }
418  for (unsigned int i = 0; i < value_n; i++) value[i] = result[i]; // cast
419  return true;
420  };
#define mycerr
bool NearestNeighbor(const coord_t *coord, value_t *value) const
Nearest-neighbor (unchecked)
enum QweakSimFieldMap::debug_t fDebug
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
EInterpolationMethod fInterpolationMethod
Interpolation method.
#define myendl
bool Linear(const coord_t *coord, value_t *value) const
Linear interpolation (unchecked)

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
unsigned int QweakSimFieldMap< value_t, value_n >::GetWrapCoordinate ( const unsigned int  dim) const
inline

Get wrapping coordinate.

Definition at line 175 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fWrap.

176  { return fWrap.at(dim); }
std::vector< size_t > fWrap
Wrap around this coordinate.
template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::InBounds ( const coord_t coord) const
inline

Return true if the coordinate is in bounds.

Definition at line 241 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Check().

241  {
242  return Check(coord);
243  };
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
unsigned int QweakSimFieldMap< value_t, value_n >::Index ( const coord_t coord) const
inline

Return the linearized index based on the point coordinates (unchecked)

Definition at line 473 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Cell(), and QweakSimFieldMap< value_t, value_n >::fNDim.

Referenced by QweakSimFieldMap< value_t, value_n >::Get(), and QweakSimFieldMap< value_t, value_n >::Set().

473  {
474  unsigned int* cell_index = new unsigned int[fNDim];
475  Cell(coord, cell_index);
476  unsigned int index = Index(cell_index);
477  delete[] cell_index;
478  return index;
479  };
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
unsigned int QweakSimFieldMap< value_t, value_n >::Index ( const unsigned int *  cell_index) const
inline

Return the linearized index based on the cell indices (unchecked)

Return the linearized index based on the cell indices (unchecked)

Parameters
cell_indexIndex in each dimension
Returns
Linearized index

Definition at line 683 of file QweakSimFieldMap.hh.

685 {
686  unsigned int linear_index = 0;
687  // Loop over dimensions
688  for (unsigned int dim = 0; dim < fNDim; dim++)
689  // ... and use the stored extents for every dimensions
690  linear_index += fExtent[dim] * cell_index[dim];
691  return linear_index;
692 }
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t , unsigned int value_n>
unsigned int QweakSimFieldMap< value_t, value_n >::Index ( const unsigned int *  cell_index,
const unsigned int  pattern 
) const
inline

Return the linearized index based on the cell indices and offset (unchecked)

Return the linearized index based on the cell indices and offset (unchecked)

Parameters
cell_indexIndex in each dimension
patternBit pattern with offsets in each dimension
Returns
Linearized index

Definition at line 701 of file QweakSimFieldMap.hh.

704 {
705  unsigned int linear_index = 0;
706  // Loop over dimensions
707  for (unsigned int dim = 0; dim < fNDim; dim++) {
708  // Add offset of one based on the bit at position dim
709  unsigned int offset = (pattern >> dim) & 0x1;
710  linear_index += fExtent[dim] * (cell_index[dim] + offset);
711  }
712  return linear_index;
713 }
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::Linear ( const coord_t coord,
value_t *  value 
) const
inlineprivate

Linear interpolation (unchecked)

Perform the multi-dimensional linear interpolation (unchecked)

Parameters
coordPoint coordinates
valueInterpolated value (reference)

Definition at line 561 of file QweakSimFieldMap.hh.

References mycout, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::GetValue().

564 {
565  // Get cell and local normalized coordinates
566  unsigned int* cell_index = new unsigned int[fNDim];
567  double* cell_local = new double[fNDim];
568  Cell(coord, cell_index, cell_local);
569  // Initialize to zero
570  for (unsigned int i = 0; i < value_n; i++) value[i] = 0.0;
571  // Calculate the interpolated value
572  // by summing the 2^fNDim = 1 << fNDim neighbor points (1U is unsigned one)
573  for (unsigned int offset = 0; offset < (1U << fNDim); offset++) {
574  value_t neighbor[value_n];
575  if (! Get(Index(cell_index,offset), neighbor)) {
576  if (fDebug >= kDebug) {
577  mycout << "QweakSimFieldMap::Linear(coord_t*,value_t*): "
578  << "neighbor " << offset << " could not be found " << myendl;
579  }
580  return false;
581  }
582  // ... with appropriate weighting factors (1 - cell_local or cell_local)
583  double fac = 1.0;
584  for (unsigned int dim = 0; dim < fNDim; dim++)
585  fac *= ((offset >> dim) & 0x1)? cell_local[dim] : (1 - cell_local[dim]);
586  // ... for all vector components
587  for (unsigned int i = 0; i < value_n; i++)
588  value[i] += fac * neighbor[i];
589  }
590  delete[] cell_index;
591  delete[] cell_local;
592  return true;
593 }
value_t Get(const unsigned int *cell_index) const
Get a single value by cell index (unchecked)
enum QweakSimFieldMap::debug_t fDebug
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::Nearest ( const coord_t coord,
unsigned int *  cell_index 
) const
inlineprivate

Return the cell index closest to the coordinate (could be above) (unchecked)

Definition at line 504 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Cell(), and QweakSimFieldMap< value_t, value_n >::fNDim.

Referenced by QweakSimFieldMap< value_t, value_n >::Set().

504  {
505  double* cell_local = new double[fNDim];
506  Cell(coord, cell_index, cell_local);
507  // Loop over all dimensions and add one if larger than 0.5
508  for (unsigned int dim = 0; dim < fNDim; dim++)
509  if (cell_local[dim] > 0.5) cell_index[dim]++;
510  delete[] cell_local;
511  };
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::NearestNeighbor ( const coord_t coord,
value_t *  value 
) const
inlineprivate

Nearest-neighbor (unchecked)

Perform the nearest-neighbor interpolation (unchecked)

Parameters
coordPoint coordinates
valueInterpolated value (reference)

Definition at line 601 of file QweakSimFieldMap.hh.

Referenced by QweakSimFieldMap< value_t, value_n >::GetValue().

604 {
605  // Get nearest cell
606  unsigned int* cell_index = new unsigned int[fNDim];
607  Nearest(coord, cell_index);
608  bool status = Get(Index(cell_index), value);
609  delete[] cell_index;
610  return status;
611 }
value_t Get(const unsigned int *cell_index) const
Get a single value by cell index (unchecked)
unsigned int fNDim
Number of dimensions in coordinates.
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)
void Nearest(const coord_t *coord, unsigned int *cell_index) const
Return the cell index closest to the coordinate (could be above) (unchecked)

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::PrintCoverage ( const unsigned int  dim)
inline

Print coverage map for all bins in one dimension.

Definition at line 209 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Cell(), QweakSimFieldMap< value_t, value_n >::Coord(), QweakSimFieldMap< value_t, value_n >::fMaximumEntries, QweakSimFieldMap< value_t, value_n >::fNDim, QweakSimFieldMap< value_t, value_n >::fSize, QweakSimFieldMap< value_t, value_n >::Get(), mycout, and myendl.

209  {
210  // Initialize coverage
211  unsigned int* cover = new unsigned int[fSize[dim]];
212  unsigned int* total = new unsigned int[fSize[dim]];
213  for (unsigned int i = 0; i < fSize[dim]; i++) {
214  cover[i] = 0;
215  total[i] = 0;
216  }
217  // Loop over all entries
218  value_t value[value_n];
219  unsigned int* cell_index = new unsigned int[fNDim];
220  for (unsigned int linear_index = 0; linear_index < fMaximumEntries; linear_index++) {
221  Cell(linear_index,cell_index);
222  total[cell_index[dim]]++;
223  if (Get(linear_index,value) && value[0] != 0)
224  cover[cell_index[dim]]++; // non-zero field
225  }
226  // Print coverage
227  coord_t* coord = new coord_t[fNDim];
228  for (size_t i = 0; i < fSize[dim]; i++) {
229  cell_index[dim] = i;
230  Coord(cell_index,coord);
231  mycout << "bin " << i << ", coord " << coord[dim] << ": "
232  << double(cover[i]) / double(total[i]) * 100 << "%"<< myendl;
233  }
234  delete[] cell_index;
235  delete[] coord;
236  delete[] cover;
237  delete[] total;
238  }
value_t Get(const unsigned int *cell_index) const
Get a single value by cell index (unchecked)
double coord_t
void Coord(const unsigned int *cell_index, coord_t *coord) const
Return the coordinates based on the cell index (unchecked)
unsigned int fMaximumEntries
Maximum number of values.
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
std::vector< size_t > fSize
Number of points in each dimension.

+ Here is the call graph for this function:

template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::ReadBinaryFile ( const std::string &  filename)
inline

Read the grid values from binary file.

Read the grid values from binary file (should be 64-bit safe, untested)

Parameters
filenameFile name
Returns
True if read successfully

Definition at line 967 of file QweakSimFieldMap.hh.

References mycerr, mycout, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap().

968 {
969  std::ifstream file(filename.c_str(), std::ios::binary);
970  if (! file.is_open()) {
971  mycerr << "QweakSimFieldMap::ReadBinaryFile(string&): "
972  << "unable to open file " << filename << "." << myendl;
973  return false;
974  }
975  // Informational message
976  mycout << "Reading binary file: ";
977  // Go to end and store length (could also use std::ios::ate)
978  file.seekg(0, std::ios::end);
979  // Go to begin and start reading template parameters
980  file.seekg(0, std::ios::beg);
981  uint32_t n, size;
982  file.read(reinterpret_cast<char*>(&n), sizeof(n));
983  file.read(reinterpret_cast<char*>(&size), sizeof(size));
984  if (n != value_n) {
985  mycerr << "QweakSimFieldMap::ReadBinaryFile(string&): "
986  << "incompatible field dimension in stream." << myendl;
987  return false;
988  }
989  if (size != sizeof(value_t)) {
990  mycerr << "QweakSimFieldMap::ReadBinaryFile(string&): "
991  << "incompatible precision of field values.." << myendl;
992  return false;
993  }
994  // Read grid parameters
995  uint32_t ndim;
996  file.read(reinterpret_cast<char*>(&ndim), sizeof(ndim));
997  SetDimensions(ndim);
998  // Read sizes
999  for (unsigned int dim = 0; dim < fNDim; dim++) {
1000  file.read(reinterpret_cast<char*>(&fMin[dim]),sizeof(fMin[dim]));
1001  file.read(reinterpret_cast<char*>(&fMax[dim]),sizeof(fMax[dim]));
1002  file.read(reinterpret_cast<char*>(&fStep[dim]),sizeof(fStep[dim]));
1003  }
1004  SetMinimumMaximumStep(fMin, fMax, fStep); // calculates fMaximumEntries
1005  // Read values
1006  uint32_t maximum_entries;
1007  file.read(reinterpret_cast<char*>(&maximum_entries),sizeof(maximum_entries));
1008  if (maximum_entries != fMaximumEntries) {
1009  mycerr << "QweakSimFieldMap::ReadBinaryFile(string&): "
1010  << "expected number of entries." << myendl;
1011  return false;
1012  }
1013  int value_size = sizeof(value_t);
1014  unsigned int entries = fValues[0].size();
1015  for (unsigned int index = 0; index < entries; index++) {
1016  // Read values
1017  for (unsigned int i = 0; i < value_n; i++) {
1018  file.read((char*)(&fValues[i][index]),value_size);
1019  }
1020  // Progress bar
1021  if (index % (entries / 10) == 0)
1022  mycout << index / (entries / 100) << "%" << std::flush;
1023  if (index % (entries / 10) != 0
1024  && index % (entries / 40) == 0)
1025  mycout << "." << std::flush;
1026  }
1027  mycout << myendl;
1028  file.close();
1029  return true;
1030 }
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycerr
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
unsigned int fMaximumEntries
Maximum number of values.
std::vector< coord_t > fStep
Step size in each dimension.
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< coord_t > fMin
Minimum in each dimension.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::ReadText ( std::istream &  stream)
inline

Read the grid from text.

Read the grid values from a text stream

Parameters
streamInput stream
Returns
True if successfully read all values

Definition at line 861 of file QweakSimFieldMap.hh.

References mycerr, mycout, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::ReadTextFile().

862 {
863  // Informational message
864  mycout << "Reading text stream: ";
865  // Read the dimensions
866  unsigned int n;
867  stream >> fNDim >> n;
868  if (n != value_n) {
869  mycerr << "QweakSimFieldMap::ReadText(istream&): "
870  << "incompatible field dimension in stream." << myendl;
871  return false;
872  }
874  // Read the grid parameters
875  for (unsigned int dim = 0; dim < fNDim; dim++)
876  stream >> dim >> fMin[dim] >> fMax[dim] >> fStep[dim];
878  // Read the grid values
879  unsigned int entries;
880  stream >> entries;
881  for (unsigned int index = 0; index < entries; index++) {
882  // Read values
883  for (unsigned int i = 0; i < value_n; i++) {
884  stream >> fValues[i][index];
885  }
886  // Progress bar
887  if (index % (entries / 10) == 0)
888  mycout << index / (entries / 100) << "%" << std::flush;
889  if (index % (entries / 10) != 0
890  && index % (entries / 40) == 0)
891  mycout << "." << std::flush;
892  }
893  // Check for end of file
894  std::string end;
895  stream >> end;
896  // Informational message
897  mycout << myendl;
898  if (end == "end") return true;
899  else {
900  mycerr << "QweakSimFieldMap::ReadText(istream&): "
901  << "stream did not end with keyword 'end'." << myendl;
902  return false;
903  }
904 }
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycerr
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
std::vector< coord_t > fStep
Step size in each dimension.
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< coord_t > fMin
Minimum in each dimension.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::ReadTextFile ( const std::string &  filename)
inline

Read the grid from text file.

Definition at line 448 of file QweakSimFieldMap.hh.

References mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::ReadText().

448  {
449  std::ifstream file(filename.c_str());
450  if (! file.is_open()) {
451  mycerr << "QweakSimFieldMap::ReadTextFile(string&): "
452  << "unable to open file " << filename << "." << myendl;
453  return false;
454  }
455  if (! ReadText(file)) {
456  mycerr << "QweakSimFieldMap::ReadTextFile(string&): "
457  << "unable to read text stream." << myendl;
458  return false;
459  }
460  file.close();
461  return true;
462  };
#define mycerr
bool ReadText(std::istream &stream)
Read the grid from text.
#define myendl

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const coord_t coord,
const value_t &  value 
)
inline

Set a single value at a coordinate (false if not possible)

Definition at line 249 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fNDim, mycerr, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::Set().

249  {
250  if (value_n != 1) {
251  mycerr << "QweakSimFieldMap::Set(coord_t&,value_t&): "
252  << "only valid for one-dimensional values." << myendl;
253  return false;
254  }
255  if (fNDim != 1) {
256  mycerr << "QweakSimFieldMap::Set(coord_t&,value_t&): "
257  << "only valid for one-dimensional grids." << myendl;
258  return false;
259  }
260  return Set(&coord, &value);
261  };
#define mycerr
#define myendl
unsigned int fNDim
Number of dimensions in coordinates.
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const coord_t coord,
const value_t &  value 
)
inline

Set a single value at a coordinate (false if not possible)

Definition at line 263 of file QweakSimFieldMap.hh.

References mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::Set().

263  {
264  if (value_n != 1) {
265  mycerr << "QweakSimFieldMap::Set(coord_t*,value_t&): "
266  << "only valid for one-dimensional values." << myendl;
267  return false;
268  }
269  return Set(coord, &value);
270  };
#define mycerr
#define myendl
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const coord_t coord,
const value_t *  value 
)
inline

Set a set of values at a coordinate (false if not possible)

Definition at line 272 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Check(), QweakSimFieldMap< value_t, value_n >::fNDim, QweakSimFieldMap< value_t, value_n >::fSize, QweakSimFieldMap< value_t, value_n >::fWrap, QweakSimFieldMap< value_t, value_n >::Index(), QweakSimFieldMap< value_t, value_n >::Nearest(), and QweakSimFieldMap< value_t, value_n >::Set().

272  {
273  unsigned int* cell_index = new unsigned int[fNDim];
274  Nearest(coord, cell_index); // nearest cell
275  if (! Check(cell_index)) return false; // out of bounds
276  bool status = true;
277  bool written = false;
278  unsigned int linear_index;
279  for (unsigned int dim = 0; dim < fNDim; dim++) {
280  // skip dimensions that are not wrapped around
281  if (fWrap[dim] == 0) continue;
282  // FIXME only one wrapping coordinate correctly supported in Set()
283  if ((cell_index[dim] < fWrap[dim]) ||
284  (fSize[dim] - cell_index[dim] - 1 < fWrap[dim])) {
285  // there are equivalent grid points
286  for (size_t wrap = 0; wrap < fWrap[dim]; wrap++) {
287  // at the minimum
288  cell_index[dim] = wrap;
289  linear_index = Index(cell_index);
290  status &= Set(linear_index, value);
291  // at the maximum
292  cell_index[dim] = fSize[dim] - wrap - 1;
293  linear_index = Index(cell_index);
294  status &= Set(linear_index, value);
295  // set flag
296  written = true;
297  }
298  }
299  }
300  if (not written) {
301  // this is an unambiguous grid point
302  linear_index = Index(cell_index);
303  status &= Set(linear_index, value);
304  }
305  delete[] cell_index;
306  return status;
307  };
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fWrap
Wrap around this coordinate.
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
std::vector< size_t > fSize
Number of points in each dimension.
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)
void Nearest(const coord_t *coord, unsigned int *cell_index) const
Return the cell index closest to the coordinate (could be above) (unchecked)

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const unsigned int  linear_index,
const value_t &  value 
)
inline

Set a single value at a linearized index (false if not possible)

Definition at line 310 of file QweakSimFieldMap.hh.

References mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::Set().

310  {
311  if (value_n != 1) {
312  mycerr << "QweakSimFieldMap::Set(unsigned int,value_t&): "
313  << "only valid for one-dimensional values." << myendl;
314  return false;
315  }
316  return Set(linear_index, &value);
317  };
#define mycerr
#define myendl
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const unsigned int  linear_index,
const value_t *  value 
)
inline

Set a set of values at a linearized index (false if not possible)

Definition at line 319 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Check(), QweakSimFieldMap< value_t, value_n >::fCurrentEntries, and QweakSimFieldMap< value_t, value_n >::fValues.

319  {
320  if (! Check(linear_index)) return false; // out of bounds
321  for (unsigned int i = 0; i < value_n; i++)
322  fValues[i][linear_index] = value[i];
323  fCurrentEntries++;
324  return true;
325  };
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
unsigned int fCurrentEntries
Number of values read in.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const unsigned int *  cell_index,
const value_t &  value 
)
inline

Set a single value at a grid point (false if out of bounds)

Definition at line 328 of file QweakSimFieldMap.hh.

References mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::Set().

328  {
329  if (value_n != 1) {
330  mycerr << "QweakSimFieldMap::Set(unsigned int*,value_t&): "
331  << "only valid for one-dimensional values." << myendl;
332  return false;
333  }
334  return Set(cell_index, &value);
335  };
#define mycerr
#define myendl
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::Set ( const unsigned int *  cell_index,
const value_t *  value 
)
inline

Set a set of values at a grid point (false if out of bounds)

Definition at line 337 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::Check(), QweakSimFieldMap< value_t, value_n >::Index(), and QweakSimFieldMap< value_t, value_n >::Set().

337  {
338  if (! Check(cell_index)) return false; // out of bounds
339  return Set(Index(cell_index),value);
340  };
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetDataReductionFactor ( const unsigned int  dim,
const unsigned int  redux 
)
inline

Set data reduction factor.

Definition at line 189 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fRedux.

Referenced by QweakSimFieldMap< value_t, value_n >::SetDataReductionFactor().

190  { fRedux.at(dim) = redux; }
std::vector< size_t > fRedux
Data reduction factor.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetDataReductionFactor ( const unsigned int  redux)
inline

Definition at line 191 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fNDim, and QweakSimFieldMap< value_t, value_n >::SetDataReductionFactor().

191  {
192  for (unsigned int dim = 0; dim < fNDim; dim++)
193  SetDataReductionFactor(dim,redux);
194  }
void SetDataReductionFactor(const unsigned int dim, const unsigned int redux)
Set data reduction factor.
unsigned int fNDim
Number of dimensions in coordinates.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetDataReductionFactor ( const std::vector< unsigned int > &  redux)
inline

Definition at line 195 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fNDim, and QweakSimFieldMap< value_t, value_n >::fRedux.

195  {
196  if (redux.size() != fNDim) return;
197  fRedux = redux;
198  }
std::vector< size_t > fRedux
Data reduction factor.
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetDebugLevel ( debug_t  debug)
inline

Set debug level.

Definition at line 81 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fDebug.

81 { fDebug = debug; };
enum QweakSimFieldMap::debug_t fDebug
template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetDimensions ( const unsigned int  ndim)
inline

Set the number of coordinate dimensions and resize vectors.

Definition at line 120 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fExtent, QweakSimFieldMap< value_t, value_n >::fMax, QweakSimFieldMap< value_t, value_n >::fMin, QweakSimFieldMap< value_t, value_n >::fNDim, QweakSimFieldMap< value_t, value_n >::fSize, QweakSimFieldMap< value_t, value_n >::fStep, QweakSimFieldMap< value_t, value_n >::fWrap, mycerr, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap(), and QweakSimFieldMap< value_t, value_n >::SetMinimumMaximumStep().

120  {
121  if (ndim == 0) {
122  mycerr << "QweakSimFieldMap::SetDimensions: "
123  << "Dimensions of the grid should be strictly positive!" << myendl;
124  return;
125  }
126  fNDim = ndim;
127  fMin.resize(fNDim); fMax.resize(fNDim); fStep.resize(fNDim); fWrap.resize(fNDim);
128  fSize.resize(fNDim); fExtent.resize(fNDim+1);
129  };
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycerr
std::vector< size_t > fExtent
std::vector< coord_t > fStep
Step size in each dimension.
#define myendl
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fWrap
Wrap around this coordinate.
std::vector< size_t > fSize
Number of points in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetInterpolationMethod ( const EInterpolationMethod  method)
inline

Set the interpolation method.

Definition at line 201 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fInterpolationMethod.

Referenced by QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap().

202  { fInterpolationMethod = method; };
EInterpolationMethod fInterpolationMethod
Interpolation method.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetMinimumMaximumStep ( const coord_t  min,
const coord_t  max,
const coord_t  step 
)
inline

Set minimum, maximum, and step size to single values.

Definition at line 131 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fNDim.

Referenced by QweakSimFieldMap< value_t, value_n >::QweakSimFieldMap().

131  {
132  SetMinimumMaximumStep(std::vector<coord_t>(fNDim,min),
133  std::vector<coord_t>(fNDim,max),
134  std::vector<coord_t>(fNDim,step));
135  };
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
unsigned int fNDim
Number of dimensions in coordinates.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetMinimumMaximumStep ( const std::vector< coord_t > &  min,
const std::vector< coord_t > &  max,
const std::vector< coord_t > &  step 
)
inline

Set minimum, maximum, and step size to different values.

Definition at line 137 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fCurrentEntries, QweakSimFieldMap< value_t, value_n >::fExtent, QweakSimFieldMap< value_t, value_n >::fMax, QweakSimFieldMap< value_t, value_n >::fMaximumEntries, QweakSimFieldMap< value_t, value_n >::fMin, QweakSimFieldMap< value_t, value_n >::fNDim, QweakSimFieldMap< value_t, value_n >::fSize, QweakSimFieldMap< value_t, value_n >::fStep, QweakSimFieldMap< value_t, value_n >::fValues, mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::SetDimensions().

139  {
140  // Set the dimensionality
141  if (min.size() != fNDim) SetDimensions(min.size());
142  // Check the dimensionality and assign boundaries and step size vectors
143  if (min.size() != fNDim || min.size() != fNDim || step.size() != fNDim) return;
144  fMin = min; fMax = max; fStep = step;
145  fExtent[0] = 1;
146  for (size_t i = 0; i < fMin.size(); i++) {
147  coord_t int_part; // safer to use modf than a direct static cast
148  coord_t frac_part = modf((fMax[i] - fMin[i]) / fStep[i] + 1.0, &int_part);
149  fSize[i] = static_cast<unsigned int>(int_part) + (frac_part > 0.5 ? 1 : 0);
150  fExtent[i+1] = fExtent[i] * fSize[i];
151  }
152  // Try resizing to allocate memory and initialize with zeroes
153  fMaximumEntries = fExtent[fNDim]; fCurrentEntries = 0; // no entries read yet
154  for (unsigned int i = 0; i < value_n; i++) {
155  try {
156  fValues[i].resize(fMaximumEntries,0);
157  } catch (std::bad_alloc) {
158  mycerr << "QweakSimFieldMap:SetMinimumMaximumStep: "
159  << "could not allocate memory for values!" << myendl;
160  }
161  }
162  };
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycerr
double coord_t
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
std::vector< size_t > fExtent
unsigned int fMaximumEntries
Maximum number of values.
std::vector< coord_t > fStep
Step size in each dimension.
#define myendl
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fSize
Number of points in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fCurrentEntries
Number of values read in.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetWrapCoordinate ( const unsigned int  dim,
const size_t  wrap = 1 
)
inline

Set wrapping coordinate.

Definition at line 178 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fWrap.

179  { fWrap.at(dim) = wrap; }
std::vector< size_t > fWrap
Wrap around this coordinate.
template<class value_t = float, unsigned int value_n = 1>
void QweakSimFieldMap< value_t, value_n >::SetWrapCoordinate ( const std::vector< size_t > &  wrap)
inline

Definition at line 180 of file QweakSimFieldMap.hh.

References QweakSimFieldMap< value_t, value_n >::fNDim, and QweakSimFieldMap< value_t, value_n >::fWrap.

180  {
181  if (wrap.size() != fNDim) return;
182  fWrap = wrap;
183  }
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fWrap
Wrap around this coordinate.
template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::WriteBinaryFile ( const std::string &  filename) const
inline

Write the grid values to binary file.

Write the grid values to binary file (should be 64-bit safe, untested)

Integer data types can be stored differently on 32-bit and 64-bit systems. Fixed-length types uint32_t and u_int32_t are provided in stdint.h and sys/types.h, respectively. The floating point types float and double will always have a length of 32 and 64 bit, per the IEEE convention.

Parameters
filenameFile name
Returns
True if written successfully

Definition at line 918 of file QweakSimFieldMap.hh.

References mycerr, mycout, and myendl.

919 {
920  std::ofstream file(filename.c_str(), std::ios::binary);
921  if (! file.is_open()) {
922  mycerr << "QweakSimFieldMap::WriteBinaryFile(string&): "
923  << "unable to open file " << filename << "." << myendl;
924  return false;
925  }
926  // Informational message
927  mycout << "Writing binary file: ";
928  // Write template parameters
929  uint32_t n = value_n; // uint32_t has length of 32 bits on any system
930  uint32_t size = sizeof(value_t);
931  file.write(reinterpret_cast<const char*>(&n), sizeof(n));
932  file.write(reinterpret_cast<const char*>(&size), sizeof(size));
933  // Write grid parameters
934  uint32_t ndim = fNDim;
935  file.write(reinterpret_cast<const char*>(&ndim), sizeof(ndim));
936  for (unsigned int dim = 0; dim < fNDim; dim++) {
937  file.write(reinterpret_cast<const char*>(&fMin[dim]),sizeof(fMin[dim]));
938  file.write(reinterpret_cast<const char*>(&fMax[dim]),sizeof(fMax[dim]));
939  file.write(reinterpret_cast<const char*>(&fStep[dim]),sizeof(fStep[dim]));
940  }
941  uint32_t entries = fValues[0].size();
942  file.write(reinterpret_cast<const char*>(&entries),sizeof(entries));
943  for (unsigned int index = 0; index < entries; index++) {
944  // Write values
945  for (unsigned int i = 0; i < value_n; i++) {
946  value_t value = fValues[i][index];
947  file.write(reinterpret_cast<const char*>(&value),sizeof(value));
948  }
949  // Progress bar
950  if (index % (entries / 10) == 0)
951  mycout << index / (entries / 100) << "%" << std::flush;
952  if (index % (entries / 10) != 0
953  && index % (entries / 40) == 0)
954  mycout << "." << std::flush;
955  }
956  mycout << myendl;
957  file.close();
958  return true;
959 }
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycerr
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
std::vector< coord_t > fStep
Step size in each dimension.
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< coord_t > fMin
Minimum in each dimension.
template<class value_t , unsigned int value_n>
bool QweakSimFieldMap< value_t, value_n >::WriteText ( std::ostream &  stream) const
inline

Write the grid as text.

Write the grid values to a text stream

Parameters
streamOutput stream

Definition at line 825 of file QweakSimFieldMap.hh.

References mycout, and myendl.

Referenced by QweakSimFieldMap< value_t, value_n >::WriteTextFile(), and QweakSimFieldMap< value_t, value_n >::WriteTextScreen().

826 {
827  // Informational message
828  mycout << "Writing text stream: ";
829  // Write the dimensions
830  stream << fNDim << "\t" << value_n << std::endl;
831  // Write the grid parameters
832  for (unsigned int dim = 0; dim < fNDim; dim++)
833  stream << dim << "\t" << fMin[dim] << "\t" << fMax[dim] << "\t" << fStep[dim] << std::endl;
834  // Write the values
835  stream << fValues[0].size() << std::endl;
836  unsigned int entries = fValues[0].size();
837  for (unsigned int index = 0; index < entries; index++) {
838  // Write values
839  for (unsigned int i = 0; i < value_n; i++) {
840  stream << fValues[i][index] << "\t";
841  }
842  stream << std::endl;
843  // Progress bar
844  if (index % (entries / 10) == 0)
845  mycout << index / (entries / 100) << "%" << std::flush;
846  if (index % (entries / 10) != 0
847  && index % (entries / 40) == 0)
848  mycout << "." << std::flush;
849  }
850  stream << "end" << std::endl;
851  mycout << myendl;
852  return true;
853 }
std::vector< coord_t > fMax
Maximum in each dimension.
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
std::vector< coord_t > fStep
Step size in each dimension.
#define myendl
#define mycout
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< coord_t > fMin
Minimum in each dimension.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::WriteTextFile ( const std::string &  filename) const
inline

Write the grid to text file.

Definition at line 429 of file QweakSimFieldMap.hh.

References mycerr, myendl, and QweakSimFieldMap< value_t, value_n >::WriteText().

429  {
430  std::ofstream file(filename.c_str());
431  if (! file.is_open()) {
432  mycerr << "QweakSimFieldMap::WriteTextFile(string&): "
433  << "unable to open file " << filename << "." << myendl;
434  return false;
435  }
436  WriteText(file);
437  file.close();
438  return true;
439  };
#define mycerr
#define myendl
bool WriteText(std::ostream &stream) const
Write the grid as text.

+ Here is the call graph for this function:

template<class value_t = float, unsigned int value_n = 1>
bool QweakSimFieldMap< value_t, value_n >::WriteTextScreen ( ) const
inline

Write the grid to screen.

Definition at line 441 of file QweakSimFieldMap.hh.

References mycout, and QweakSimFieldMap< value_t, value_n >::WriteText().

441  {
442  WriteText(mycout);
443  return true;
444  };
#define mycout
bool WriteText(std::ostream &stream) const
Write the grid as text.

+ Here is the call graph for this function:

Field Documentation

template<class value_t = float, unsigned int value_n = 1>
unsigned int QweakSimFieldMap< value_t, value_n >::fCurrentEntries
private
template<class value_t = float, unsigned int value_n = 1>
enum QweakSimFieldMap::debug_t QweakSimFieldMap< value_t, value_n >::fDebug
template<class value_t = float, unsigned int value_n = 1>
std::vector<size_t> QweakSimFieldMap< value_t, value_n >::fExtent
private

Linear extent between neighbor points in each dimension (e.g. for the least significant index this will be 1, for the next index the number of points in the first index, etc...)

Definition at line 104 of file QweakSimFieldMap.hh.

Referenced by QweakSimFieldMap< value_t, value_n >::SetDimensions(), and QweakSimFieldMap< value_t, value_n >::SetMinimumMaximumStep().

template<class value_t = float, unsigned int value_n = 1>
EInterpolationMethod QweakSimFieldMap< value_t, value_n >::fInterpolationMethod
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<coord_t> QweakSimFieldMap< value_t, value_n >::fMax
private
template<class value_t = float, unsigned int value_n = 1>
unsigned int QweakSimFieldMap< value_t, value_n >::fMaximumEntries
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<coord_t> QweakSimFieldMap< value_t, value_n >::fMin
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<size_t> QweakSimFieldMap< value_t, value_n >::fRedux
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<size_t> QweakSimFieldMap< value_t, value_n >::fSize
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<coord_t> QweakSimFieldMap< value_t, value_n >::fStep
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<value_t> QweakSimFieldMap< value_t, value_n >::fValues[value_n]
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<size_t> QweakSimFieldMap< value_t, value_n >::fWrap
private

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