QwAnalysis
QwInterpolator< value_t, value_n > Class Template Reference

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

#include <QwInterpolator.h>

+ Inheritance diagram for QwInterpolator< value_t, value_n >:

Public Member Functions

 QwInterpolator (const unsigned int ndim=1)
 Constructor with number of dimensions. More...
 
 QwInterpolator (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...
 
 QwInterpolator (const std::string &filename)
 Constructor with file name. More...
 
virtual ~QwInterpolator ()
 Destructor. 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 EQwInterpolationMethod method)
 Set the interpolation method. More...
 
EQwInterpolationMethod 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...
 

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...
 
unsigned int * f_cell_index
 Pre-allocated cell index. More...
 
double * f_cell_local
 Pre-allocated local coordinates. More...
 
EQwInterpolationMethod fInterpolationMethod
 Interpolation method. More...
 

Detailed Description

template<class value_t = float, unsigned int value_n = 1>
class QwInterpolator< 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 QwInterpolator<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 QwInterpolator.h.

Constructor & Destructor Documentation

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

Constructor with number of dimensions.

Definition at line 56 of file QwInterpolator.h.

56  {
57  SetDimensions(ndim);
59  f_cell_index = new unsigned int[fNDim];
60  f_cell_local = new double[fNDim];
61  };
unsigned int fNDim
Number of dimensions in coordinates.
void SetInterpolationMethod(const EQwInterpolationMethod method)
Set the interpolation method.
unsigned int * f_cell_index
Pre-allocated cell index.
double * f_cell_local
Pre-allocated local coordinates.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
template<class value_t = float, unsigned int value_n = 1>
QwInterpolator< value_t, value_n >::QwInterpolator ( 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 63 of file QwInterpolator.h.

65  {
66  SetDimensions(min.size());
68  SetMinimumMaximumStep(min,max,step);
69  f_cell_index = new unsigned int[fNDim];
70  f_cell_local = new double[fNDim];
71  };
unsigned int fNDim
Number of dimensions in coordinates.
void SetInterpolationMethod(const EQwInterpolationMethod method)
Set the interpolation method.
unsigned int * f_cell_index
Pre-allocated cell index.
double * f_cell_local
Pre-allocated local coordinates.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
template<class value_t = float, unsigned int value_n = 1>
QwInterpolator< value_t, value_n >::QwInterpolator ( const std::string &  filename)
inline

Constructor with file name.

Definition at line 73 of file QwInterpolator.h.

73  {
74  ReadBinaryFile(filename);
76  f_cell_index = new unsigned int[fNDim];
77  f_cell_local = new double[fNDim];
78  };
unsigned int fNDim
Number of dimensions in coordinates.
void SetInterpolationMethod(const EQwInterpolationMethod method)
Set the interpolation method.
bool ReadBinaryFile(const std::string &filename)
Read the grid values from binary file.
unsigned int * f_cell_index
Pre-allocated cell index.
double * f_cell_local
Pre-allocated local coordinates.
template<class value_t = float, unsigned int value_n = 1>
virtual QwInterpolator< value_t, value_n >::~QwInterpolator ( )
inlinevirtual

Destructor.

Definition at line 80 of file QwInterpolator.h.

80  {
81  delete[] f_cell_index;
82  delete[] f_cell_local;
83  };
unsigned int * f_cell_index
Pre-allocated cell index.
double * f_cell_local
Pre-allocated local coordinates.

Member Function Documentation

template<class value_t , unsigned int value_n>
void QwInterpolator< 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 609 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::Index(), QwInterpolator< float, 4 >::Nearest(), and QwInterpolator< float, 4 >::PrintCoverage().

614 {
615  // Normalize coordinate (positive, with integers on grid points)
616  double norm_coord = (coord - fMin[dim]) / fStep[dim];
617  // Split in fractional and integer part
618  double frac_part, int_part;
619  frac_part = modf(norm_coord, &int_part);
620  cell_local = frac_part;
621  cell_index = static_cast<int>(int_part); // cast to integer
622  // Wrap index
623  if (fWrap[dim] > 0) cell_index %= (fSize[dim] - fWrap[dim]);
624 }
std::vector< coord_t > fStep
Step size in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.
std::vector< size_t > fSize
Number of points in each dimension.
std::vector< size_t > fWrap
Wrap around this coordinate.

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
void QwInterpolator< 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 633 of file QwInterpolator.h.

637 {
638  // Get cell index and local coordinates in each dimension
639  for (unsigned int dim = 0; dim < fNDim; dim++)
640  Cell(coord[dim], cell_index[dim], cell_local[dim], dim);
641 }
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 QwInterpolator< 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 649 of file QwInterpolator.h.

652 {
653  // Get cell index and ignore local coordinates
654  Cell(coord, cell_index, f_cell_local);
655 }
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)
double * f_cell_local
Pre-allocated local coordinates.
template<class value_t , unsigned int value_n>
void QwInterpolator< 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 663 of file QwInterpolator.h.

666 {
667  // This does not work with unsigned int, because that is always >= 0 and wraps
668  for (int dim = fNDim-1; dim >= 0; dim--) {
669  cell_index[dim] = linear_index / fExtent[dim];
670  linear_index -= cell_index[dim] * fExtent[dim];
671  }
672 }
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t , unsigned int value_n>
bool QwInterpolator< 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 525 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::GetValue(), QwInterpolator< float, 4 >::InBounds(), and QwInterpolator< float, 4 >::Set().

526 {
527  for (unsigned int dim = 0; dim < fNDim; dim++)
528  if (fWrap[dim] == 0 && (coord[dim] < fMin[dim] || coord[dim] > fMax[dim]))
529  return false;
530  // Otherwise
531  return true;
532 }
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< coord_t > fMax
Maximum in each dimension.
std::vector< size_t > fWrap
Wrap around this coordinate.

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
bool QwInterpolator< 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 540 of file QwInterpolator.h.

541 {
542  for (unsigned int dim = 0; dim < fNDim; dim++)
543  if (cell_index[dim] >= fSize[dim])
544  return false;
545  // Otherwise
546  return true;
547 }
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 QwInterpolator< 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 555 of file QwInterpolator.h.

556 {
557  if (linear_index >= fMaximumEntries)
558  return false;
559  // Otherwise
560  return true;
561 }
unsigned int fMaximumEntries
Maximum number of values.
template<class value_t , unsigned int value_n>
void QwInterpolator< 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 680 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::PrintCoverage(), and QwMatrixLookup::WriteTrajMatrix().

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

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
void QwInterpolator< 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 694 of file QwInterpolator.h.

697 {
698  Cell(linear_index,f_cell_index);
699  Coord(f_cell_index,coord);
700 }
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 * f_cell_index
Pre-allocated cell index.
void Coord(const unsigned int *cell_index, coord_t *coord) const
Return the coordinates based on the cell index (unchecked)
template<class value_t = float, unsigned int value_n = 1>
value_t QwInterpolator< value_t, value_n >::Get ( const unsigned int *  cell_index) const
inlineprivate

Get a single value by cell index (unchecked)

Definition at line 454 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::PrintCoverage().

454  {
455  if (value_n != 1) return 0; // only for one-dimensional values
456  return fValues[0][Index(cell_index)];
457  };
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.

+ Here is the caller graph for this function:

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

Get a single value by linearized index (unchecked)

Definition at line 460 of file QwInterpolator.h.

460  {
461  return fValues[0][index];
462  };
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 464 of file QwInterpolator.h.

464  {
465  for (unsigned int i = 0; i < value_n; i++)
466  value[i] = fValues[i][index];
467  return true;
468  };
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 QwInterpolator< value_t, value_n >::GetCurrentEntries ( ) const
inline

Get the current number of entries.

Definition at line 178 of file QwInterpolator.h.

Referenced by QwMatrixLookup::LoadTrajMatrix(), QwMagneticField::ReadFieldMapStream(), and QwMatrixLookup::WriteTrajMatrix().

178 { return fCurrentEntries; };
unsigned int fCurrentEntries
Number of values read in.

+ Here is the caller graph for this function:

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

Get data reduction factor.

Definition at line 192 of file QwInterpolator.h.

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

Get the interpolation method.

Definition at line 210 of file QwInterpolator.h.

Referenced by QwMagneticField::GetInterpolationMethod().

211  { return fInterpolationMethod; };
EQwInterpolationMethod fInterpolationMethod
Interpolation method.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
coord_t QwInterpolator< value_t, value_n >::GetMaximum ( const unsigned int  dim) const
inline

Get maximum in dimension.

Definition at line 172 of file QwInterpolator.h.

Referenced by QwMatrixLookup::Bridge().

172 { return fMax.at(dim); };
std::vector< coord_t > fMax
Maximum in each dimension.

+ Here is the caller graph for this function:

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

Get the maximum number of entries.

Definition at line 176 of file QwInterpolator.h.

Referenced by QwMatrixLookup::LoadTrajMatrix(), QwMagneticField::ReadFieldMapStream(), and QwMatrixLookup::WriteTrajMatrix().

176 { return fMaximumEntries; };
unsigned int fMaximumEntries
Maximum number of values.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
coord_t QwInterpolator< value_t, value_n >::GetMinimum ( const unsigned int  dim) const
inline

Get minimum in dimension.

Definition at line 170 of file QwInterpolator.h.

Referenced by QwMatrixLookup::Bridge().

170 { return fMin.at(dim); };
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>
coord_t QwInterpolator< value_t, value_n >::GetStepSize ( const unsigned int  dim) const
inline

Get minimum in dimension.

Definition at line 174 of file QwInterpolator.h.

Referenced by QwMatrixLookup::Bridge().

174 { return fStep.at(dim); };
std::vector< coord_t > fStep
Step size in each dimension.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
value_t QwInterpolator< 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 329 of file QwInterpolator.h.

Referenced by QwMatrixLookup::Bridge(), QwMagneticField::GetFieldValue(), and QwInterpolator< float, 4 >::GetValue().

329  {
330  if (value_n != 1) return false; // only for one-dimensional values
331  if (fNDim != 1) return false; // only for one-dimensional grids
332  value_t value;
333  if (GetValue(&coord, &value)) return value;
334  else return 0; // interpolation failed
335  };
unsigned int fNDim
Number of dimensions in coordinates.
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
value_t QwInterpolator< 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 337 of file QwInterpolator.h.

337  {
338  if (value_n != 1) return 0; // only for one-dimensional values
339  value_t value;
340  if (GetValue(coord, &value)) return value;
341  else return 0; // interpolation failed
342  };
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 344 of file QwInterpolator.h.

344  {
345  if (value_n != 1) return false; // only for one-dimensional values
346  return GetValue(coord, &value);
347  };
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 349 of file QwInterpolator.h.

349  {
350  for (unsigned int i = 0; i < value_n; i++) value[i] = 0.0; // zero
351  if (! Check(coord)) return false; // out of bounds
352  value_t result[value_n]; // we need a local copy of type value_t
353  switch (fInterpolationMethod) {
354  // return interpolation status
355  case kMultiLinear:
356  if (! Linear(coord, result)) return false;
357  break;
358  case kNearestNeighbor:
359  if (! NearestNeighbor(coord, result)) return false;
360  break;
361  default: return false;
362  }
363  for (unsigned int i = 0; i < value_n; i++) value[i] = result[i]; // cast
364  return true;
365  };
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
bool Linear(const coord_t *coord, value_t *value) const
Linear interpolation (unchecked)
bool NearestNeighbor(const coord_t *coord, value_t *value) const
Nearest-neighbor (unchecked)
EQwInterpolationMethod fInterpolationMethod
Interpolation method.
template<class value_t = float, unsigned int value_n = 1>
unsigned int QwInterpolator< value_t, value_n >::GetWrapCoordinate ( const unsigned int  dim) const
inline

Get wrapping coordinate.

Definition at line 181 of file QwInterpolator.h.

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

Return true if the coordinate is in bounds.

Definition at line 245 of file QwInterpolator.h.

Referenced by QwMatrixLookup::Bridge(), and QwMagneticField::GetFieldValue().

245  {
246  return Check(coord);
247  };
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.

+ Here is the caller graph for this function:

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

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

Definition at line 406 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::Get(), QwInterpolator< float, 4 >::Index(), and QwInterpolator< float, 4 >::Set().

406  {
407  Cell(coord, f_cell_index);
408  return Index(f_cell_index);
409  };
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 * f_cell_index
Pre-allocated cell index.
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 , unsigned int value_n>
unsigned int QwInterpolator< 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 569 of file QwInterpolator.h.

571 {
572  unsigned int linear_index = 0;
573  // Loop over dimensions
574  for (unsigned int dim = 0; dim < fNDim; dim++)
575  // ... and use the stored extents for every dimensions
576  linear_index += fExtent[dim] * cell_index[dim];
577  return linear_index;
578 }
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t , unsigned int value_n>
unsigned int QwInterpolator< 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 587 of file QwInterpolator.h.

590 {
591  unsigned int linear_index = 0;
592  // Loop over dimensions
593  for (unsigned int dim = 0; dim < fNDim; dim++) {
594  // Add offset of one based on the bit at position dim
595  unsigned int offset = (pattern >> dim) & 0x1;
596  linear_index += fExtent[dim] * (cell_index[dim] + offset);
597  }
598  return linear_index;
599 }
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
template<class value_t, unsigned int value_n>
bool QwInterpolator< 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 480 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::GetValue().

483 {
484  // Get cell and local normalized coordinates
485  Cell(coord, f_cell_index, f_cell_local);
486  // Initialize to zero
487  for (unsigned int i = 0; i < value_n; i++) value[i] = 0.0;
488  // Calculate the interpolated value
489  // by summing the 2^fNDim = 1 << fNDim neighbor points (1U is unsigned one)
490  for (unsigned int offset = 0; offset < (1U << fNDim); offset++) {
491  value_t neighbor[value_n];
492  if (! Get(Index(f_cell_index,offset), neighbor)) return false;
493  // ... with appropriate weighting factors (1 - cell_local or cell_local)
494  double fac = 1.0;
495  for (unsigned int dim = 0; dim < fNDim; dim++)
496  fac *= ((offset >> dim) & 0x1)? f_cell_local[dim] : (1 - f_cell_local[dim]);
497  // ... for all vector components
498  for (unsigned int i = 0; i < value_n; i++)
499  value[i] += fac * neighbor[i];
500  }
501  return true;
502 }
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 * f_cell_index
Pre-allocated cell index.
value_t Get(const unsigned int *cell_index) const
Get a single value by cell index (unchecked)
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)
double * f_cell_local
Pre-allocated local coordinates.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< 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 434 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::Set().

434  {
435  Cell(coord, cell_index, f_cell_local);
436  // Loop over all dimensions and add one if larger than 0.5
437  for (unsigned int dim = 0; dim < fNDim; dim++)
438  if (f_cell_local[dim] > 0.5) cell_index[dim]++;
439  };
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)
double * f_cell_local
Pre-allocated local coordinates.

+ Here is the caller graph for this function:

template<class value_t, unsigned int value_n>
bool QwInterpolator< 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 510 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::GetValue().

513 {
514  // Get nearest cell
515  Nearest(coord, f_cell_index);
516  return Get(Index(f_cell_index), value);
517 }
unsigned int * f_cell_index
Pre-allocated cell index.
value_t Get(const unsigned int *cell_index) const
Get a single value by cell index (unchecked)
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 QwInterpolator< value_t, value_n >::PrintCoverage ( const unsigned int  dim)
inline

Print coverage map for all bins in one dimension.

Definition at line 215 of file QwInterpolator.h.

Referenced by QwMagneticField::ReadFieldMapStream().

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

+ Here is the caller graph for this function:

template<class value_t , unsigned int value_n>
bool QwInterpolator< 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 837 of file QwInterpolator.h.

References QwLog::flush(), mycout, and myendl.

Referenced by QwInterpolator< float, 4 >::QwInterpolator(), and QwMagneticField::ReadFieldMap().

838 {
839  std::ifstream file(filename.c_str(), std::ios::binary);
840  if (! file.is_open()) return false;
841  // Informational message
842  mycout << "Reading binary file: ";
843  // Go to end and store length (could also use std::ios::ate)
844  file.seekg(0, std::ios::end);
845  // Go to begin and start reading template parameters
846  file.seekg(0, std::ios::beg);
847  uint32_t n, size;
848  file.read(reinterpret_cast<char*>(&n), sizeof(n));
849  file.read(reinterpret_cast<char*>(&size), sizeof(size));
850  if (n != value_n) return false; // not same dimensionality
851  if (size != sizeof(value_t)) return false; // not same type
852  // Read grid parameters
853  uint32_t ndim;
854  file.read(reinterpret_cast<char*>(&ndim), sizeof(ndim));
855  SetDimensions(ndim);
856  // Read sizes
857  for (unsigned int dim = 0; dim < fNDim; dim++) {
858  file.read(reinterpret_cast<char*>(&fMin[dim]),sizeof(fMin[dim]));
859  file.read(reinterpret_cast<char*>(&fMax[dim]),sizeof(fMax[dim]));
860  file.read(reinterpret_cast<char*>(&fStep[dim]),sizeof(fStep[dim]));
861  }
862  SetMinimumMaximumStep(fMin, fMax, fStep); // calculates fMaximumEntries
863  // Read values
864  uint32_t maximum_entries;
865  file.read(reinterpret_cast<char*>(&maximum_entries),sizeof(maximum_entries));
866  if (maximum_entries != fMaximumEntries) return false; // not expected number of entries
867  int value_size = sizeof(value_t);
868  unsigned int entries = fValues[0].size();
869  for (unsigned int index = 0; index < entries; index++) {
870  // Read values
871  for (unsigned int i = 0; i < value_n; i++) {
872  file.read((char*)(&fValues[i][index]),value_size);
873  }
874  // Progress bar
875  if (index % (entries / 10) == 0)
876  mycout << index / (entries / 100) << "%" << QwLog::flush;
877  if (index % (entries / 10) != 0
878  && index % (entries / 40) == 0)
879  mycout << "." << QwLog::flush;
880  }
881  mycout << myendl;
882  file.close();
883  return true;
884 }
std::vector< coord_t > fStep
Step size in each dimension.
unsigned int fMaximumEntries
Maximum number of values.
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
#define myendl
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycout
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
static std::ostream & flush(std::ostream &)
Flush the streams.
Definition: QwLog.cc:319

+ 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 QwInterpolator< 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 743 of file QwInterpolator.h.

References QwLog::flush(), mycout, and myendl.

Referenced by QwInterpolator< float, 4 >::ReadTextFile().

744 {
745  // Informational message
746  mycout << "Reading text stream: ";
747  // Read the dimensions
748  unsigned int n;
749  stream >> fNDim >> n;
750  if (n != value_n) return false; // not same dimensionality
752  // Read the grid parameters
753  for (unsigned int dim = 0; dim < fNDim; dim++)
754  stream >> dim >> fMin[dim] >> fMax[dim] >> fStep[dim];
756  // Read the grid values
757  unsigned int entries;
758  stream >> entries;
759  for (unsigned int index = 0; index < entries; index++) {
760  // Read values
761  for (unsigned int i = 0; i < value_n; i++) {
762  stream >> fValues[i][index];
763  }
764  // Progress bar
765  if (index % (entries / 10) == 0)
766  mycout << index / (entries / 100) << "%" << QwLog::flush;
767  if (index % (entries / 10) != 0
768  && index % (entries / 40) == 0)
769  mycout << "." << QwLog::flush;
770  }
771  // Check for end of file
772  std::string end;
773  stream >> end;
774  // Informational message
775  mycout << myendl;
776  if (end == "end") return true;
777  else return false;
778 }
std::vector< coord_t > fStep
Step size in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
#define myendl
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycout
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
static std::ostream & flush(std::ostream &)
Flush the streams.
Definition: QwLog.cc:319

+ 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>
bool QwInterpolator< value_t, value_n >::ReadTextFile ( const std::string &  filename)
inline

Read the grid from text file.

Definition at line 389 of file QwInterpolator.h.

389  {
390  std::ifstream file(filename.c_str());
391  if (! file.is_open()) return false;
392  if (! ReadText(file)) return false;
393  file.close();
394  return true;
395  };
bool ReadText(std::istream &stream)
Read the grid from text.
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 253 of file QwInterpolator.h.

Referenced by QwMatrixLookup::LoadTrajMatrix(), QwMagneticField::ReadFieldMapStream(), QwInterpolator< float, 4 >::Set(), and QwMatrixLookup::WriteTrajMatrix().

253  {
254  if (value_n != 1) return false; // only for one-dimensional values
255  if (fNDim != 1) return false; // only for one-dimensional grids
256  return Set(&coord, &value);
257  };
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 QwInterpolator< 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 259 of file QwInterpolator.h.

259  {
260  if (value_n != 1) return false; // only for one-dimensional values
261  return Set(coord, &value);
262  };
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 264 of file QwInterpolator.h.

264  {
265  Nearest(coord, f_cell_index); // nearest cell
266  if (! Check(f_cell_index)) return false; // out of bounds
267  bool status = true;
268  bool written = false;
269  unsigned int linear_index;
270  for (unsigned int dim = 0; dim < fNDim; dim++) {
271  // skip dimensions that are not wrapped around
272  if (fWrap[dim] == 0) continue;
273  // FIXME only one wrapping coordinate correctly supported in Set()
274  if ((f_cell_index[dim] < fWrap[dim]) ||
275  (fSize[dim] - f_cell_index[dim] - 1 < fWrap[dim])) {
276  // there are equivalent grid points
277  for (size_t wrap = 0; wrap < fWrap[dim]; wrap++) {
278  // at the minimum
279  f_cell_index[dim] = wrap;
280  linear_index = Index(f_cell_index);
281  status &= Set(linear_index, value);
282  // at the maximum
283  f_cell_index[dim] = fSize[dim] - wrap - 1;
284  linear_index = Index(f_cell_index);
285  status &= Set(linear_index, value);
286  // set flag
287  written = true;
288  }
289  }
290  }
291  if (not written) {
292  // this is an unambiguous grid point
293  linear_index = Index(f_cell_index);
294  status &= Set(linear_index, value);
295  }
296  return status;
297  };
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)
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
std::vector< size_t > fSize
Number of points in each dimension.
unsigned int * f_cell_index
Pre-allocated cell index.
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)
std::vector< size_t > fWrap
Wrap around this coordinate.
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 300 of file QwInterpolator.h.

300  {
301  if (value_n != 1) return false; // only for one-dimensional values
302  return Set(linear_index, &value);
303  };
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 305 of file QwInterpolator.h.

305  {
306  if (! Check(linear_index)) return false; // out of bounds
307  for (unsigned int i = 0; i < value_n; i++)
308  fValues[i][linear_index] = value[i];
309  fCurrentEntries++;
310  return true;
311  };
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
unsigned int fCurrentEntries
Number of values read in.
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 314 of file QwInterpolator.h.

314  {
315  if (value_n != 1) return false; // only for one-dimensional values
316  return Set(cell_index, &value);
317  };
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
template<class value_t = float, unsigned int value_n = 1>
bool QwInterpolator< 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 319 of file QwInterpolator.h.

319  {
320  if (! Check(cell_index)) return false; // out of bounds
321  return Set(Index(cell_index),value);
322  };
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)
template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< value_t, value_n >::SetDataReductionFactor ( const unsigned int  dim,
const unsigned int  redux 
)
inline

Set data reduction factor.

Definition at line 195 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::SetDataReductionFactor().

196  { 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 QwInterpolator< value_t, value_n >::SetDataReductionFactor ( const unsigned int  redux)
inline

Definition at line 197 of file QwInterpolator.h.

197  {
198  for (unsigned int dim = 0; dim < fNDim; dim++)
199  SetDataReductionFactor(dim,redux);
200  }
unsigned int fNDim
Number of dimensions in coordinates.
void SetDataReductionFactor(const unsigned int dim, const unsigned int redux)
Set data reduction factor.
template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< value_t, value_n >::SetDataReductionFactor ( const std::vector< unsigned int > &  redux)
inline

Definition at line 201 of file QwInterpolator.h.

201  {
202  if (redux.size() != fNDim) return;
203  fRedux = redux;
204  }
unsigned int fNDim
Number of dimensions in coordinates.
std::vector< size_t > fRedux
Data reduction factor.
template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< value_t, value_n >::SetDimensions ( const unsigned int  ndim)
inline

Set the number of coordinate dimensions and resize vectors.

Definition at line 128 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::QwInterpolator(), and QwInterpolator< float, 4 >::SetMinimumMaximumStep().

128  {
129  if (ndim == 0) {
130  mycerr << "Dimensions of the grid should be strictly positive!" << myendl;
131  return;
132  }
133  fNDim = ndim;
134  fMin.resize(fNDim); fMax.resize(fNDim); fStep.resize(fNDim); fWrap.resize(fNDim);
135  fSize.resize(fNDim); fExtent.resize(fNDim+1);
136  };
std::vector< coord_t > fStep
Step size in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.
std::vector< size_t > fExtent
unsigned int fNDim
Number of dimensions in coordinates.
#define myendl
std::vector< coord_t > fMax
Maximum in each dimension.
std::vector< size_t > fSize
Number of points in each dimension.
std::vector< size_t > fWrap
Wrap around this coordinate.
#define mycerr

+ Here is the caller graph for this function:

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

Set the interpolation method.

Definition at line 207 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::QwInterpolator(), QwMatrixLookup::QwMatrixLookup(), and QwMagneticField::SetInterpolationMethod().

208  { fInterpolationMethod = method; };
EQwInterpolationMethod fInterpolationMethod
Interpolation method.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< 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 138 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::QwInterpolator(), QwMatrixLookup::QwMatrixLookup(), and QwInterpolator< float, 4 >::SetMinimumMaximumStep().

138  {
139  SetMinimumMaximumStep(std::vector<coord_t>(fNDim,min),
140  std::vector<coord_t>(fNDim,max),
141  std::vector<coord_t>(fNDim,step));
142  };
unsigned int fNDim
Number of dimensions in coordinates.
static const double min
Definition: QwUnits.h:76
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< 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 144 of file QwInterpolator.h.

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

Set wrapping coordinate.

Definition at line 184 of file QwInterpolator.h.

Referenced by QwMagneticField::ReadFieldMap().

185  { fWrap.at(dim) = wrap; }
std::vector< size_t > fWrap
Wrap around this coordinate.

+ Here is the caller graph for this function:

template<class value_t = float, unsigned int value_n = 1>
void QwInterpolator< value_t, value_n >::SetWrapCoordinate ( const std::vector< size_t > &  wrap)
inline

Definition at line 186 of file QwInterpolator.h.

186  {
187  if (wrap.size() != fNDim) return;
188  fWrap = wrap;
189  }
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 QwInterpolator< 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 792 of file QwInterpolator.h.

References QwLog::flush(), mycout, and myendl.

Referenced by QwMagneticField::WriteBinaryFile().

793 {
794  std::ofstream file(filename.c_str(), std::ios::binary);
795  if (! file.is_open()) return false;
796  // Informational message
797  mycout << "Writing binary file: ";
798  // Write template parameters
799  uint32_t n = value_n; // uint32_t has length of 32 bits on any system
800  uint32_t size = sizeof(value_t);
801  file.write(reinterpret_cast<const char*>(&n), sizeof(n));
802  file.write(reinterpret_cast<const char*>(&size), sizeof(size));
803  // Write grid parameters
804  uint32_t ndim = fNDim;
805  file.write(reinterpret_cast<const char*>(&ndim), sizeof(ndim));
806  for (unsigned int dim = 0; dim < fNDim; dim++) {
807  file.write(reinterpret_cast<const char*>(&fMin[dim]),sizeof(fMin[dim]));
808  file.write(reinterpret_cast<const char*>(&fMax[dim]),sizeof(fMax[dim]));
809  file.write(reinterpret_cast<const char*>(&fStep[dim]),sizeof(fStep[dim]));
810  }
811  uint32_t entries = fValues[0].size();
812  file.write(reinterpret_cast<const char*>(&entries),sizeof(entries));
813  for (unsigned int index = 0; index < entries; index++) {
814  // Write values
815  for (unsigned int i = 0; i < value_n; i++) {
816  value_t value = fValues[i][index];
817  file.write(reinterpret_cast<const char*>(&value),sizeof(value));
818  }
819  // Progress bar
820  if (index % (entries / 10) == 0)
821  mycout << index / (entries / 100) << "%" << QwLog::flush;
822  if (index % (entries / 10) != 0
823  && index % (entries / 40) == 0)
824  mycout << "." << QwLog::flush;
825  }
826  mycout << myendl;
827  file.close();
828  return true;
829 }
std::vector< coord_t > fStep
Step size in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
#define myendl
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycout
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
static std::ostream & flush(std::ostream &)
Flush the streams.
Definition: QwLog.cc:319

+ 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 QwInterpolator< 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 707 of file QwInterpolator.h.

References QwLog::flush(), mycout, and myendl.

Referenced by QwInterpolator< float, 4 >::WriteTextFile(), and QwInterpolator< float, 4 >::WriteTextScreen().

708 {
709  // Informational message
710  mycout << "Writing text stream: ";
711  // Write the dimensions
712  stream << fNDim << "\t" << value_n << std::endl;
713  // Write the grid parameters
714  for (unsigned int dim = 0; dim < fNDim; dim++)
715  stream << dim << "\t" << fMin[dim] << "\t" << fMax[dim] << "\t" << fStep[dim] << std::endl;
716  // Write the values
717  stream << fValues[0].size() << std::endl;
718  unsigned int entries = fValues[0].size();
719  for (unsigned int index = 0; index < entries; index++) {
720  // Write values
721  for (unsigned int i = 0; i < value_n; i++) {
722  stream << fValues[i][index] << "\t";
723  }
724  stream << std::endl;
725  // Progress bar
726  if (index % (entries / 10) == 0)
727  mycout << index / (entries / 100) << "%" << QwLog::flush;
728  if (index % (entries / 10) != 0
729  && index % (entries / 40) == 0)
730  mycout << "." << QwLog::flush;
731  }
732  stream << "end" << std::endl;
733  mycout << myendl;
734  return true;
735 }
std::vector< coord_t > fStep
Step size in each dimension.
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
#define myendl
std::vector< coord_t > fMax
Maximum in each dimension.
#define mycout
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
static std::ostream & flush(std::ostream &)
Flush the streams.
Definition: QwLog.cc:319

+ 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>
bool QwInterpolator< value_t, value_n >::WriteTextFile ( const std::string &  filename) const
inline

Write the grid to text file.

Definition at line 374 of file QwInterpolator.h.

Referenced by QwMagneticField::WriteTextFile().

374  {
375  std::ofstream file(filename.c_str());
376  if (! file.is_open()) return false;
377  WriteText(file);
378  file.close();
379  return true;
380  };
bool WriteText(std::ostream &stream) const
Write the grid as text.

+ Here is the caller graph for this function:

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

Write the grid to screen.

Definition at line 382 of file QwInterpolator.h.

382  {
383  WriteText(std::cout);
384  return true;
385  };
bool WriteText(std::ostream &stream) const
Write the grid as text.

Field Documentation

template<class value_t = float, unsigned int value_n = 1>
unsigned int* QwInterpolator< value_t, value_n >::f_cell_index
private
template<class value_t = float, unsigned int value_n = 1>
double* QwInterpolator< value_t, value_n >::f_cell_local
private
template<class value_t = float, unsigned int value_n = 1>
unsigned int QwInterpolator< value_t, value_n >::fCurrentEntries
private
template<class value_t = float, unsigned int value_n = 1>
std::vector<size_t> QwInterpolator< 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 106 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::SetDimensions(), and QwInterpolator< float, 4 >::SetMinimumMaximumStep().

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

Table with pointers to arrays of values.

Definition at line 109 of file QwInterpolator.h.

Referenced by QwInterpolator< float, 4 >::Get(), QwInterpolator< float, 4 >::Set(), and QwInterpolator< float, 4 >::SetMinimumMaximumStep().

template<class value_t = float, unsigned int value_n = 1>
std::vector<size_t> QwInterpolator< value_t, value_n >::fWrap
private

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