SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SSAGES::GridBase< T > Class Template Referenceabstract

Base class for Grids. More...

#include <GridBase.h>

Inheritance diagram for SSAGES::GridBase< T >:
Inheritance graph
[legend]

Public Member Functions

void syncGrid ()
 Sync the grid.
 
size_t GetDimension () const
 Get the dimension. More...
 
const std::vector< int > GetNumPoints () const
 Get the number of points for all dimensions. More...
 
int GetNumPoints (size_t dim) const
 Get the number of points for a specific dimension. More...
 
const std::vector< double > GetLower () const
 Return the lower edges of the Grid. More...
 
double GetLower (size_t dim) const
 Get the lower edge for a specific dimension. More...
 
const std::vector< double > GetUpper () const
 Return the upper edges of the Grid. More...
 
double GetUpper (size_t dim) const
 Get the upper edge for a specific dimension. More...
 
const std::vector< bool > & GetPeriodic () const
 Return the periodicity of the Grid. More...
 
bool GetPeriodic (size_t dim) const
 Get the periodicity in a specific dimension. More...
 
size_t size () const
 Get the size of the internal storage vector. More...
 
T * data ()
 Get pointer to the internal data storage vector. More...
 
T const * data () const
 Get pointer to const of the internal data storage vector. More...
 
std::vector< int > GetIndices (const std::vector< double > &x) const
 Return the Grid indices for a given point. More...
 
int GetIndex (double x) const
 Return linear interpolation on a coordinate. More...
 
std::vector< double > GetCoordinates (const std::vector< int > &indices)
 ! Return the distance to adjacent grid center points from given coordinates More...
 
double GetCoordinate (int index)
 Return center point of 1d-grid. More...
 
const T & at (const std::vector< int > &indices) const
 Access Grid element read-only. More...
 
T & at (const std::vector< int > &indices)
 Access Grid element read/write. More...
 
template<typename R >
const T & at (std::initializer_list< R > &&x) const
 Const access of Grid element via initializer list. More...
 
template<typename R >
T & at (std::initializer_list< R > &&x)
 Access Grid element via initializer list. More...
 
const T & at (int index) const
 Access 1d Grid by index, read-only. More...
 
T & at (int index)
 Access 1d Grid by index, read-write. More...
 
const T & at (const std::vector< double > &x) const
 Access Grid element pertaining to a specific point – read-only. More...
 
T & at (const std::vector< double > &x)
 Access Grid element pertaining to a specific point – read/write. More...
 
const T & at (double x) const
 Access 1d-Grid by point - read-only. More...
 
T & at (double x)
 Access 1d-Grid by point - read-write. More...
 
const T & operator[] (const std::vector< int > &indices) const
 Access Grid element per [] read-only. More...
 
T & operator[] (const std::vector< int > &indices)
 Access Grid element per [] read-write. More...
 
template<typename R >
const T & operator[] (std::initializer_list< R > &&x) const
 Const access of Grid element via initializer list. More...
 
template<typename R >
T & operator[] (std::initializer_list< R > &&x)
 Access Grid element via initializer list. More...
 
const T & operator[] (int index) const
 Access 1d-Grid per [] operator, read-only. More...
 
T & operator[] (int index)
 Access 1d-Grid per [] operator, read-write. More...
 
const T & operator[] (const std::vector< double > &x) const
 Access Grid element pertaining to a specific point per [] read-only. More...
 
T & operator[] (const std::vector< double > &x)
 Access Grid element pertaining to a specific point per [] read-write. More...
 
const T & operator[] (double x) const
 Access 1d-Grid via specific point, read-only. More...
 
T & operator[] (double x)
 Access 1d-Grid via specific point, read-write. More...
 
virtual ~GridBase ()
 Destructor.
 

Protected Member Functions

std::vector< int > wrapIndices (const std::vector< int > &indices) const
 Wrap the index around periodic boundaries. More...
 
virtual size_t mapTo1d (const std::vector< int > &indices) const =0
 This function needs to be implemented by child classes. More...
 
 GridBase (std::vector< int > numPoints, std::vector< double > lower, std::vector< double > upper, std::vector< bool > isPeriodic)
 Constructor. More...
 

Protected Attributes

std::vector< T > data_
 Internal storage of the data.
 
size_t dimension_
 Dimension of the grid.
 
std::vector< int > numPoints_
 Number of points in each dimension.
 
std::pair< std::vector< double >, std::vector< double > > edges_
 Edges of the Grid in each dimension.
 
std::vector< bool > isPeriodic_
 Periodicity of the Grid.
 

Detailed Description

template<typename T>
class SSAGES::GridBase< T >

Base class for Grids.

Template Parameters
Typeof data stored on the grid

Base class for all grids. Currently, these are 'Grid' and 'Histogram'.

Definition at line 38 of file GridBase.h.

Constructor & Destructor Documentation

◆ GridBase()

template<typename T >
SSAGES::GridBase< T >::GridBase ( std::vector< int >  numPoints,
std::vector< double >  lower,
std::vector< double >  upper,
std::vector< bool >  isPeriodic 
)
inlineprotected

Constructor.

Parameters
numPointsNumber of grid points in each dimension.
lowerLower edges of the grid.
upperUpper edges of the grid.
isPeriodicBools specifying the periodicity in the respective dimension.

The constructor is protected by design. This makes sure that only child classes of GridBase are constructed.

Definition at line 104 of file GridBase.h.

108  : dimension_(numPoints.size()),
109  numPoints_(numPoints),
110  edges_(std::pair< std::vector<double>, std::vector<double> >(lower, upper)),
111  isPeriodic_(isPeriodic)
112  {
113  // Check that vector sizes are correct
114  if (edges_.first.size() != dimension_ ||
115  edges_.second.size() != dimension_) {
116  throw std::invalid_argument("Size of vector containing upper or "
117  "lower edges, does not match size of vector containing "
118  "number of grid points.");
119  }
120  if (isPeriodic_.size() == 0) {
121  // Default: Non-periodic in all dimensions
122  isPeriodic.resize(dimension_, false);
123  } else if (isPeriodic_.size() != dimension_) {
124  throw std::invalid_argument("Size of vector isPeriodic does not "
125  "match size of vector containing number of grid points.");
126  }
127  for(size_t i=0 ; i < isPeriodic_.size() ; i++)
128  {
129  if(isPeriodic_[i])
130  {
131  if(numPoints_[i] <= 1)
132  {
133  throw std::invalid_argument("A periodic grid is incompatible with a grid size of 1.");
134  }
135  double spacing = (edges_.second[i] - edges_.first[i]) / (numPoints_[i]-1);
136  edges_.first[i] -= spacing/2;
137  edges_.second[i] += spacing/2;
138  }
139  }
140 
141  }
std::vector< bool > isPeriodic_
Periodicity of the Grid.
Definition: GridBase.h:54
std::pair< std::vector< double >, std::vector< double > > edges_
Edges of the Grid in each dimension.
Definition: GridBase.h:51
std::vector< int > numPoints_
Number of points in each dimension.
Definition: GridBase.h:48
size_t dimension_
Dimension of the grid.
Definition: GridBase.h:45

References SSAGES::GridBase< T >::dimension_, SSAGES::GridBase< T >::edges_, SSAGES::GridBase< T >::isPeriodic_, and SSAGES::GridBase< T >::numPoints_.

Member Function Documentation

◆ at() [1/10]

template<typename T >
T& SSAGES::GridBase< T >::at ( const std::vector< double > &  x)
inline

Access Grid element pertaining to a specific point – read/write.

Parameters
xVector of doubles specifying a point.
Returns
Reference to the value at the given coordinates.

This function is provided for convenience. It is identical to GridBase::at(GridBase::GetIndices(x)).

Definition at line 650 of file GridBase.h.

651  {
652  return at(GetIndices(x));
653  }
std::vector< int > GetIndices(const std::vector< double > &x) const
Return the Grid indices for a given point.
Definition: GridBase.h:372
const T & at(const std::vector< int > &indices) const
Access Grid element read-only.
Definition: GridBase.h:546

References SSAGES::GridBase< T >::at(), and SSAGES::GridBase< T >::GetIndices().

Here is the call graph for this function:

◆ at() [2/10]

template<typename T >
const T& SSAGES::GridBase< T >::at ( const std::vector< double > &  x) const
inline

Access Grid element pertaining to a specific point – read-only.

Parameters
xVector of doubles specifying a point.
Returns
Const reference of the value at the given coordinates.

This function is provided for convenience. It is identical to GridBase::at(GridBase::GetIndices(x)).

Definition at line 637 of file GridBase.h.

638  {
639  return at(GetIndices(x));
640  }

References SSAGES::GridBase< T >::at(), and SSAGES::GridBase< T >::GetIndices().

Here is the call graph for this function:

◆ at() [3/10]

template<typename T >
T& SSAGES::GridBase< T >::at ( const std::vector< int > &  indices)
inline

Access Grid element read/write.

Parameters
indicesVector of integers specifying the grid point.
Returns
Reference to value at the specified grid point.

Definition at line 562 of file GridBase.h.

563  {
564  return const_cast<T&>(static_cast<const GridBase<T>* >(this)->at(indices));
565  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ at() [4/10]

template<typename T >
const T& SSAGES::GridBase< T >::at ( const std::vector< int > &  indices) const
inline

Access Grid element read-only.

Parameters
indicesVector of integers specifying the grid point.
Returns
const reference of the value stored at the given grid point.

In non-periodic dimensions, the index needs to be in the interval [-1, numPoints]. Grid::at(-1) accessed the underflow bin, Grid::at(numPoints) accesses the overflow bin.

In periodic dimensions, the index may take any integer value and will be mapped back to the interval [0, numPoints-1]. Thus, Grid::at(-1) will access the same value as Grid::at(numPoints-1).

Definition at line 546 of file GridBase.h.

547  {
548  // Check that indices are in bound.
549  if (indices.size() != GetDimension()) {
550  throw std::invalid_argument("Dimension of indices does not match "
551  "dimension of the grid.");
552  }
553 
554  return data_.at(mapTo1d(indices));
555  }
virtual size_t mapTo1d(const std::vector< int > &indices) const =0
This function needs to be implemented by child classes.
size_t GetDimension() const
Get the dimension.
Definition: GridBase.h:195
std::vector< T > data_
Internal storage of the data.
Definition: GridBase.h:42

References SSAGES::GridBase< T >::data_, SSAGES::GridBase< T >::GetDimension(), and SSAGES::GridBase< T >::mapTo1d().

Referenced by SSAGES::GridBase< T >::at(), SSAGES::GridBase< T >::GetCoordinate(), SSAGES::GridBase< T >::GetIndex(), SSAGES::GridBase< T >::operator[](), SSAGES::ANN::PostIntegration(), and SSAGES::GridBase< T >::syncGrid().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ at() [5/10]

template<typename T >
T& SSAGES::GridBase< T >::at ( double  x)
inline

Access 1d-Grid by point - read-write.

Parameters
xAccess grid point pertaining to this value.
Returns
Reference to the value pertaining to the specified coordinate.
Note
This function can only be used for 1d Grids.

Definition at line 679 of file GridBase.h.

680  {
681  return const_cast<T&>(static_cast<const GridBase<T>* >(this)->at(x));
682  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ at() [6/10]

template<typename T >
const T& SSAGES::GridBase< T >::at ( double  x) const
inline

Access 1d-Grid by point - read-only.

Parameters
xAccess grid point pertaining to this value.
Returns
Const reference to the value pertaining to the specified coordinate.
Note
This function can only be used for 1d-Grids.

Definition at line 663 of file GridBase.h.

664  {
665  if (dimension_ != 1) {
666  throw std::invalid_argument("Only 1d-Grids can be accessed with a "
667  "single float as the specified point.");
668  }
669  return at({x});
670  }

References SSAGES::GridBase< T >::at(), and SSAGES::GridBase< T >::dimension_.

Here is the call graph for this function:

◆ at() [7/10]

template<typename T >
T& SSAGES::GridBase< T >::at ( int  index)
inline

Access 1d Grid by index, read-write.

Parameters
indexIndex specifying the grid point.
Returns
Reference of value at the given grid point.
Note
This function can only be used for 1d-Grids.

Definition at line 624 of file GridBase.h.

625  {
626  return const_cast<T&>(static_cast<const GridBase<T>* >(this)->at(index));
627  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ at() [8/10]

template<typename T >
const T& SSAGES::GridBase< T >::at ( int  index) const
inline

Access 1d Grid by index, read-only.

Parameters
indexIndex specifying the grid point.
Returns
Const reference of value at the given index.
Note
This function can only be used for 1d-Grids.

Definition at line 608 of file GridBase.h.

609  {
610  if (dimension_ != 1) {
611  throw std::invalid_argument("Only 1d-Grids can be accessed with a "
612  "single integer as the index.");
613  }
614  return at({index});
615  }

References SSAGES::GridBase< T >::at(), and SSAGES::GridBase< T >::dimension_.

Here is the call graph for this function:

◆ at() [9/10]

template<typename T >
template<typename R >
T& SSAGES::GridBase< T >::at ( std::initializer_list< R > &&  x)
inline

Access Grid element via initializer list.

Template Parameters
RDatatype in the initializer list
Parameters
xinitializer list
Returns
Reference to value at the specified point.

This function avoids abiguity if at() is called with a brace-enclosed initializer list. The template parameter makes sure that this function can be called with either ints, specifying a grid point, or doubles, specifying coordinates in space, inside the initializer list.

Definition at line 596 of file GridBase.h.

597  {
598  return at(static_cast<std::vector<R> >(x));
599  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ at() [10/10]

template<typename T >
template<typename R >
const T& SSAGES::GridBase< T >::at ( std::initializer_list< R > &&  x) const
inline

Const access of Grid element via initializer list.

Template Parameters
RDatatype in the initializer list
Parameters
xinitializer list
Returns
Const reference to value at the specified point.

This function avoids abiguity if at() is called with a brace-enclosed initializer list. The template parameter makes sure that this function can be called with either ints, specifying a grid point, or doubles, specifying coordinates in space, inside the initializer list.

Definition at line 579 of file GridBase.h.

580  {
581  return at(static_cast<std::vector<R> >(x));
582  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ data() [1/2]

template<typename T >
T* SSAGES::GridBase< T >::data ( )
inline

Get pointer to the internal data storage vector.

Returns
Pointer to data in the internal storage vector.

It is discouraged to directly access the internal data storage. It might, however be necessary. For example when communicating the data over MPI.

Definition at line 338 of file GridBase.h.

339  {
340  return data_.data();
341  }

References SSAGES::GridBase< T >::data_.

Referenced by SSAGES::BFS::ProjectBias(), SSAGES::ANN::ReadBias(), SSAGES::ANN::TrainNetwork(), and SSAGES::ANN::WriteBias().

Here is the caller graph for this function:

◆ data() [2/2]

template<typename T >
T const* SSAGES::GridBase< T >::data ( ) const
inline

Get pointer to const of the internal data storage vector.

Returns
Const pointer to data in the internal storage vector.

It is discouraged to directly access the internal data storage. It might, however be necessary. For example when communicating data over MPI.

Definition at line 350 of file GridBase.h.

351  {
352  return data_.data();
353  }

References SSAGES::GridBase< T >::data_.

◆ GetCoordinate()

template<typename T >
double SSAGES::GridBase< T >::GetCoordinate ( int  index)
inline

Return center point of 1d-grid.

Parameters
indexIndex of the 1d grid.
Returns
Coordinate in real/CV space.
Note
This function is only available for 1d grids.

Definition at line 528 of file GridBase.h.

529  {
530  return GetCoordinates({index}).at(0);
531  }
std::vector< double > GetCoordinates(const std::vector< int > &indices)
! Return the distance to adjacent grid center points from given coordinates
Definition: GridBase.h:503

References SSAGES::GridBase< T >::at(), and SSAGES::GridBase< T >::GetCoordinates().

Here is the call graph for this function:

◆ GetCoordinates()

template<typename T >
std::vector<double> SSAGES::GridBase< T >::GetCoordinates ( const std::vector< int > &  indices)
inline

! Return the distance to adjacent grid center points from given coordinates

Return coordinates of the grid center points

Parameters
indicesGrid indices specifying a grid point.
Returns
Vector of double specifying the position of the grid point.

The grid is a discretization of real or cv space. Thus, each grid point is associated with an interval of the underlying space. This function returns the center point of this interval.

Definition at line 503 of file GridBase.h.

504  {
505  if (indices.size() != dimension_) {
506  throw std::invalid_argument(
507  "Grid indices specified for GetCoordinates() do not have "
508  "the same dimensionality as the grid.");
509  }
510 
511  std::vector<double> v(dimension_);
512 
513  for (size_t i = 0; i < dimension_; ++i) {
514  double spacing = (edges_.second[i] - edges_.first[i]) / (numPoints_[i]);
515  v.at(i) = edges_.first[i] + (indices[i] + 0.5)*spacing;
516  }
517 
518  return v;
519  }

References SSAGES::GridBase< T >::dimension_, SSAGES::GridBase< T >::edges_, and SSAGES::GridBase< T >::numPoints_.

Referenced by SSAGES::GridBase< T >::GetCoordinate().

Here is the caller graph for this function:

◆ GetDimension()

template<typename T >
size_t SSAGES::GridBase< T >::GetDimension ( ) const
inline

◆ GetIndex()

template<typename T >
int SSAGES::GridBase< T >::GetIndex ( double  x) const
inline

Return linear interpolation on a coordinate.

   \param x Point in space.
   \return Linearly interpolated value of grid at that point.
    This function performs a n-linear interpolation on the grid.
    The formula is bilinear/trilinear interpolation generalized
    to N-dimensional space. 
 double GetInterpolated(const std::vector<double> &x)
{
    double interpvalue = 0;
    int tempindex;
    for (size_t i = 0 ; i < pow(2,dimension_) ; ++i)
    {

        std::vector<double> shiftedvector = x;
        tempindex = i;

        double accumulatedweight = 1;
        for(size_t j = 0; j < dimension_ ; ++j)
        {
            double spacing = (edges_.second[j] - edges_.first[j]) / numPoints_[j];
            shiftedvector[j] += ((tempindex%2)-0.5)*spacing;
            tempindex = tempindex/2;
        }

        std::vector<int> shiftedindices = GetIndices(shiftedvector);
        std::vector<double> shiftedcenters = GetCoordinates(shiftedindices);

        for(size_t j = 0; j < dimension_ ; ++j)
        {   
            double spacing = (edges_.second[j] - edges_.first[j]) / numPoints_[j];

Handle Edges if(shiftedindices[j] == -1) { accumulatedweight *= ((std::abs(x[j]-GetCoordinates(GetIndices(x))[j])/spacing)); shiftedvector[j] += spacing/2; } else if (shiftedindices[j] == numPoints_[j]) { accumulatedweight *= ((std::abs(x[j]-GetCoordinates(GetIndices(x))[j])/spacing)); shiftedvector[j] -= spacing/2; } else { Handle Periodicity if(std::abs(x[j]-shiftedcenters[j]) > spacing) accumulatedweight *= (1-(std::abs(std::abs(x[j]-shiftedcenters[j]) - (edges_.second[j] - edges_.first[j]))/spacing)); else accumulatedweight *= (1-(std::abs(x[j]-shiftedcenters[j])/spacing)); } }
interpvalue += accumulatedweight*at(shiftedvector); } return interpvalue; } Return the Grid index for a one-dimensional grid.

Parameters
xPoint in space.
Returns
Grid index to which the point pertains.

Return the Grid index pertaining to the given point in space. This function is for convenience when accessing 1d-Grids. For higher-dimensional grids, x needs to be a vector of doubles.

Definition at line 472 of file GridBase.h.

473  {
474  if (dimension_ != 1) {
475  throw std::invalid_argument("1d Grid index can only be accessed for "
476  "1d-Grids can be accessed with a.");
477  }
478  return GetIndices({x}).at(0);
479  }

References SSAGES::GridBase< T >::at(), SSAGES::GridBase< T >::dimension_, and SSAGES::GridBase< T >::GetIndices().

Here is the call graph for this function:

◆ GetIndices()

template<typename T >
std::vector<int> SSAGES::GridBase< T >::GetIndices ( const std::vector< double > &  x) const
inline

Return the Grid indices for a given point.

Parameters
xPoint in space.
Returns
Indices of the grid point to which the point in space pertains.

The grid discretizes the continuous space. For a given point in this continuous space, this function will return the indices of the grid point covering the point in space.

If the grid is non-periodic in a given dimension and x is lower than the lower edge in this dimension, the function will return -1, the index of the underflow bin. Similarly, it will return numPoints, the index of the overflow bin, if x is larger than the upper edge.

In periodic dimensions, the index can take any integer value and will be wrapped to the interval [0, numPoints).

Definition at line 372 of file GridBase.h.

373  {
374  // Check that input vector has the correct dimensionality
375  if (x.size() != dimension_) {
376  throw std::invalid_argument("Specified point has a larger "
377  "dimensionality than the grid.");
378  }
379 
380  std::vector<int> indices(dimension_);
381  for (size_t i = 0; i < dimension_; ++i)
382  {
383  double xpos = x.at(i);
384  if (!GetPeriodic(i))
385  {
386  if (xpos < edges_.first[i])
387  {
388  indices.at(i) = -1;
389  continue;
390  }
391  else if (xpos > edges_.second[i])
392  {
393  indices.at(i) = numPoints_[i];
394  continue;
395  }
396 
397  }
398 
399  // To make sure, the value is rounded in the correct direction.
400  double spacing = (edges_.second[i] - edges_.first[i]) / (numPoints_[i]);
401  indices.at(i) = std::floor( (xpos - edges_.first[i]) / spacing);
402  }
403 
404  return wrapIndices(indices);
405  }
const std::vector< bool > & GetPeriodic() const
Return the periodicity of the Grid.
Definition: GridBase.h:295
std::vector< int > wrapIndices(const std::vector< int > &indices) const
Wrap the index around periodic boundaries.
Definition: GridBase.h:62

References SSAGES::GridBase< T >::dimension_, SSAGES::GridBase< T >::edges_, SSAGES::GridBase< T >::GetPeriodic(), SSAGES::GridBase< T >::numPoints_, and SSAGES::GridBase< T >::wrapIndices().

Referenced by SSAGES::GridBase< T >::at(), and SSAGES::GridBase< T >::GetIndex().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetLower() [1/2]

template<typename T >
const std::vector<double> SSAGES::GridBase< T >::GetLower ( ) const
inline

Return the lower edges of the Grid.

Returns
Vector containing the lower edges of the grid.

Definition at line 231 of file GridBase.h.

232  {
233  std::vector<double> lower(dimension_);
234  for(size_t i = 0; i<dimension_;++i)
235  if(GetPeriodic(i)) {lower[i] = edges_.first[i] - ((edges_.first[i] - edges_.second[i]) / (numPoints_[i]))/2;}
236  else {lower[i] = edges_.first[i];}
237  return lower;
238  }

References SSAGES::GridBase< T >::dimension_, SSAGES::GridBase< T >::edges_, SSAGES::GridBase< T >::GetPeriodic(), and SSAGES::GridBase< T >::numPoints_.

Referenced by SSAGES::ANN::PostIntegration(), SSAGES::BasisEvaluator::UpdateBias(), and SSAGES::Grid< T >::WriteToFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetLower() [2/2]

template<typename T >
double SSAGES::GridBase< T >::GetLower ( size_t  dim) const
inline

Get the lower edge for a specific dimension.

Parameters
dimIndex of the dimension.
Returns
Value of the lower edge in the requested dimension.
Note
The first dimension has the index 0.

Definition at line 247 of file GridBase.h.

248  {
249  if (dim >= GetDimension()) {
250  std::cerr << "Warning! Lower edge requested for a dimension larger "
251  "than the grid dimensionality!\n";
252  return 0.0;
253  }
254  if(GetPeriodic(dim)){return edges_.first[dim] - ((edges_.first[dim] - edges_.second[dim]) / (numPoints_[dim]))/2;}
255  else{return edges_.first[dim];}
256  }

References SSAGES::GridBase< T >::edges_, SSAGES::GridBase< T >::GetDimension(), SSAGES::GridBase< T >::GetPeriodic(), and SSAGES::GridBase< T >::numPoints_.

Here is the call graph for this function:

◆ GetNumPoints() [1/2]

template<typename T >
const std::vector<int> SSAGES::GridBase< T >::GetNumPoints ( ) const
inline

Get the number of points for all dimensions.

Returns
Vector of ints containing the number of grid points for each dimension.

Definition at line 205 of file GridBase.h.

206  {
207  return numPoints_;
208  }

References SSAGES::GridBase< T >::numPoints_.

Referenced by SSAGES::BFS::Build(), SSAGES::BasisEvaluator::UpdateBias(), SSAGES::BasisEvaluator::UpdateCoeff(), and SSAGES::Grid< T >::WriteToFile().

Here is the caller graph for this function:

◆ GetNumPoints() [2/2]

template<typename T >
int SSAGES::GridBase< T >::GetNumPoints ( size_t  dim) const
inline

Get the number of points for a specific dimension.

Parameters
dimIndex of the dimension.
Returns
Number of grid points in the requested dimension.
Note
The first dimension uses the index 0.

Definition at line 217 of file GridBase.h.

218  {
219  if (dim >= GetDimension()) {
220  std::cerr << "Warning! Grid size requested for a dimension larger "
221  "than the grid dimensionality!\n";
222  return 0;
223  }
224  return numPoints_.at(dim);
225  }

References SSAGES::GridBase< T >::GetDimension(), and SSAGES::GridBase< T >::numPoints_.

Here is the call graph for this function:

◆ GetPeriodic() [1/2]

template<typename T >
const std::vector<bool>& SSAGES::GridBase< T >::GetPeriodic ( ) const
inline

Return the periodicity of the Grid.

Returns
Vector of bools. The values are True (False ) if the grid is periodic (non-periodic) in the given dimension.

Definition at line 295 of file GridBase.h.

296  {
297  return isPeriodic_;
298  }

References SSAGES::GridBase< T >::isPeriodic_.

Referenced by SSAGES::GridBase< T >::GetIndices(), SSAGES::GridBase< T >::GetLower(), SSAGES::GridBase< T >::GetPeriodic(), SSAGES::GridBase< T >::GetUpper(), and SSAGES::GridBase< T >::wrapIndices().

Here is the caller graph for this function:

◆ GetPeriodic() [2/2]

template<typename T >
bool SSAGES::GridBase< T >::GetPeriodic ( size_t  dim) const
inline

Get the periodicity in a specific dimension.

Parameters
dimIndex of the dimension.
Returns
True (False ) if the grid is periodic (non-periodic) in the specified dimension.
Note
The dimensions are indexed starting with 0.

Definition at line 308 of file GridBase.h.

309  {
310  if (dim >= GetDimension()) {
311  std::cerr << "Warning! Periodicity requested for a dimension larger "
312  "than the grid dimensionality!\n";
313  return false;
314  }
315  return GetPeriodic().at(dim);
316  }

References SSAGES::GridBase< T >::GetDimension(), and SSAGES::GridBase< T >::GetPeriodic().

Here is the call graph for this function:

◆ GetUpper() [1/2]

template<typename T >
const std::vector<double> SSAGES::GridBase< T >::GetUpper ( ) const
inline

Return the upper edges of the Grid.

Returns
Vector containing the upper edges of the grid.

Definition at line 262 of file GridBase.h.

263  {
264  std::vector<double> upper(dimension_);
265  for(size_t i = 0; i<dimension_;++i)
266  if(GetPeriodic(i)) {upper[i] = edges_.second[i] + ((edges_.first[i] - edges_.second[i]) / (numPoints_[i]))/2;}
267  else {upper[i] = edges_.second[i];}
268  return upper;
269  }

References SSAGES::GridBase< T >::dimension_, SSAGES::GridBase< T >::edges_, SSAGES::GridBase< T >::GetPeriodic(), and SSAGES::GridBase< T >::numPoints_.

Referenced by SSAGES::ANN::PostIntegration(), SSAGES::BasisEvaluator::UpdateBias(), and SSAGES::Grid< T >::WriteToFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetUpper() [2/2]

template<typename T >
double SSAGES::GridBase< T >::GetUpper ( size_t  dim) const
inline

Get the upper edge for a specific dimension.

Parameters
dimIndex of the dimension.
Returns
Value of the upper edge in the given dimension.
Note
The dimensions are indexed starting with 0.

Definition at line 279 of file GridBase.h.

280  {
281  if (dim >= GetDimension()) {
282  std::cerr << "Warning! Upper edge requested for a dimension larger "
283  "than the grid dimensionality!\n";
284  return 0.0;
285  }
286  if(GetPeriodic(dim)){return edges_.second[dim] + ((edges_.first[dim] - edges_.second[dim]) / (numPoints_[dim]))/2;}
287  else{return edges_.second[dim];}
288  }

References SSAGES::GridBase< T >::edges_, SSAGES::GridBase< T >::GetDimension(), SSAGES::GridBase< T >::GetPeriodic(), and SSAGES::GridBase< T >::numPoints_.

Here is the call graph for this function:

◆ mapTo1d()

template<typename T >
virtual size_t SSAGES::GridBase< T >::mapTo1d ( const std::vector< int > &  indices) const
protectedpure virtual

This function needs to be implemented by child classes.

Parameters
indicesThe indices specifying the grid point.
Returns
Index of the grid point in the 1d storage vector.

mapTo1d maps the indices onto a single index to access the data stored in the data_ vector. This mapping will be different for different grid implementations.

Implemented in SSAGES::Grid< T >, SSAGES::Grid< double >, SSAGES::Grid< Vector >, SSAGES::Grid< int >, SSAGES::Grid< unsigned int >, SSAGES::Grid< std::vector< double > >, and SSAGES::Grid< Eigen::VectorXd >.

Referenced by SSAGES::GridBase< T >::at().

Here is the caller graph for this function:

◆ operator[]() [1/10]

template<typename T >
T& SSAGES::GridBase< T >::operator[] ( const std::vector< double > &  x)
inline

Access Grid element pertaining to a specific point per [] read-write.

Parameters
xVector of doubles specifying the point in space.
Returns
Reference to the value pertaining to the given coordinates.

Definition at line 777 of file GridBase.h.

778  {
779  return at(x);
780  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [2/10]

template<typename T >
const T& SSAGES::GridBase< T >::operator[] ( const std::vector< double > &  x) const
inline

Access Grid element pertaining to a specific point per [] read-only.

Parameters
xVector of doubles specifying the point in space.
Returns
Const reference to the value pertaining to the given coordinates.

Definition at line 767 of file GridBase.h.

768  {
769  return at(x);
770  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [3/10]

template<typename T >
T& SSAGES::GridBase< T >::operator[] ( const std::vector< int > &  indices)
inline

Access Grid element per [] read-write.

Parameters
indicesVector of integers specifying the grid point.
Returns
Reference to value at the specified grid point.

Definition at line 699 of file GridBase.h.

700  {
701  return at(indices);
702  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [4/10]

template<typename T >
const T& SSAGES::GridBase< T >::operator[] ( const std::vector< int > &  indices) const
inline

Access Grid element per [] read-only.

Parameters
indicesVector of integers specifying the grid point.
Returns
Const reference to value to the given grid point.

Definition at line 689 of file GridBase.h.

690  {
691  return at(indices);
692  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [5/10]

template<typename T >
T& SSAGES::GridBase< T >::operator[] ( double  x)
inline

Access 1d-Grid via specific point, read-write.

Parameters
xPoint specifying the desired Grid point.
Returns
Reference to value at the specified coordinate.
Note
This function can only be used for 1d grids.

Definition at line 801 of file GridBase.h.

802  {
803  return at(x);
804  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [6/10]

template<typename T >
const T& SSAGES::GridBase< T >::operator[] ( double  x) const
inline

Access 1d-Grid via specific point, read-only.

Parameters
xPoint specifying the desired Grid point.
Returns
Const reference to the value at the specified coordinate.
Note
This function can only be used for 1d grids.

Definition at line 789 of file GridBase.h.

790  {
791  return at(x);
792  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [7/10]

template<typename T >
T& SSAGES::GridBase< T >::operator[] ( int  index)
inline

Access 1d-Grid per [] operator, read-write.

Parameters
indexIndex of the grid point.
Returns
Reference to the value at the given grid point.
Note
This function can only be used for 1d grids.

Definition at line 757 of file GridBase.h.

758  {
759  return at(index);
760  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [8/10]

template<typename T >
const T& SSAGES::GridBase< T >::operator[] ( int  index) const
inline

Access 1d-Grid per [] operator, read-only.

Parameters
indexIndex of the grid point.
Returns
Const reference to the value at the given grid point.
Note
This function can only be used for 1d grids.

Definition at line 745 of file GridBase.h.

746  {
747  return at(index);
748  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [9/10]

template<typename T >
template<typename R >
T& SSAGES::GridBase< T >::operator[] ( std::initializer_list< R > &&  x)
inline

Access Grid element via initializer list.

Template Parameters
RDatatype in the initializer list
Parameters
xinitializer list
Returns
Reference to the value at the specified point.

This function avoids abiguity if operator[] is called with a brace-enclosed initializer list.

Example: grid[{0,1}]

Definition at line 733 of file GridBase.h.

734  {
735  return at(static_cast<std::vector<R> >(x));
736  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ operator[]() [10/10]

template<typename T >
template<typename R >
const T& SSAGES::GridBase< T >::operator[] ( std::initializer_list< R > &&  x) const
inline

Const access of Grid element via initializer list.

Template Parameters
RDatatype in the initializer list
Parameters
xinitializer list
Returns
Const reference to the value at the specified point.

This function avoids abiguity if operator[] is called with a brace-enclosed initializer list.

Example: grid[{0,1}] or grid[{-1.23, 4.2, 0.0}]

Definition at line 716 of file GridBase.h.

717  {
718  return at(static_cast<std::vector<R> >(x));
719  }

References SSAGES::GridBase< T >::at().

Here is the call graph for this function:

◆ size()

template<typename T >
size_t SSAGES::GridBase< T >::size ( ) const
inline

Get the size of the internal storage vector.

Returns
Size of the internal storage vector.

This function returns the size of the internal storage vector. This is also the total number of grid points including the over/underflow bins in case of a histogram.

Definition at line 326 of file GridBase.h.

327  {
328  return data_.size();
329  }

References SSAGES::GridBase< T >::data_.

Referenced by SSAGES::ANN::ANN(), SSAGES::Grid< T >::BuildGrid(), SSAGES::CFF::CFF(), and SSAGES::ANN::TrainNetwork().

Here is the caller graph for this function:

◆ wrapIndices()

template<typename T >
std::vector<int> SSAGES::GridBase< T >::wrapIndices ( const std::vector< int > &  indices) const
inlineprotected

Wrap the index around periodic boundaries.

Parameters
indicesCurrent indices of grid to wrap
Returns
Resulting indices after wrapping

Definition at line 62 of file GridBase.h.

63  {
64  std::vector<int> newIndices(indices);
65  for (size_t i = 0; i < dimension_; ++i) {
66  if (!GetPeriodic(i)) {
67  continue;
68  }
69 
70  int index = indices.at(i);
71  while (index < 0) { index += numPoints_[i]; }
72  while (index >= numPoints_[i]) { index -= numPoints_[i]; }
73  if(index == numPoints_[i]-1)
74  index = 0;
75  newIndices.at(i) = index;
76  }
77 
78  return newIndices;
79  }

References SSAGES::GridBase< T >::dimension_, SSAGES::GridBase< T >::GetPeriodic(), and SSAGES::GridBase< T >::numPoints_.

Referenced by SSAGES::GridBase< T >::GetIndices().

Here is the call graph for this function:
Here is the caller graph for this function:

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