VAPoR  3.0.0
Public Member Functions | List of all members
VAPoR::LayeredGrid Class Reference

#include <LayeredGrid.h>

Inheritance diagram for VAPoR::LayeredGrid:
VAPoR::RegularGrid

Public Member Functions

 LayeredGrid (const size_t bs[3], const size_t min[3], const size_t max[3], const double extents[6], const bool periodic[3], float **blks, float **coords, int varying_dim)
 
 LayeredGrid (const size_t bs[3], const size_t min[3], const size_t max[3], const double extents[6], const bool periodic[3], float **blks, float **coords, int varying_dim, float missing_value)
 
virtual ~LayeredGrid ()
 
float GetValue (double x, double y, double z) const
 
virtual void GetUserExtents (double extents[6]) const
 
virtual void GetBoundingBox (const size_t min[3], const size_t max[3], double extents[6]) const
 
virtual void GetEnclosingRegion (const double minu[3], const double maxu[3], size_t min[3], size_t max[3]) const
 
int GetUserCoordinates (size_t i, size_t j, size_t k, double *x, double *y, double *z) const
 
void GetIJKIndex (double x, double y, double z, size_t *i, size_t *j, size_t *k) const
 
void GetIJKIndexFloor (double x, double y, double z, size_t *i, size_t *j, size_t *k) const
 
int Reshape (const size_t min[3], const size_t max[3], const bool periodic[3])
 
bool InsideGrid (double x, double y, double z) const
 
virtual void SetPeriodic (const bool periodic[3])
 
virtual void GetMinCellExtents (double *x, double *y, double *z) const
 
float ** GetCoordBlks () const
 
- Public Member Functions inherited from VAPoR::RegularGrid
 RegularGrid (const size_t bs[3], const size_t min[3], const size_t max[3], const double extents[6], const bool periodic[3], float **blks)
 
 RegularGrid (const size_t bs[3], const size_t min[3], const size_t max[3], const double extents[6], const bool periodic[3], float **blks, float missing_value)
 
virtual ~RegularGrid ()
 
float & AccessIJK (size_t i, size_t j, size_t k) const
 
virtual void GetBoundingBox (const std::vector< size_t > &min, const std::vector< size_t > &max, std::vector< double > &minu, std::vector< double > &maxu) const
 
virtual void GetIJKMinMax (size_t min[3], size_t max[3]) const
 
virtual void GetIJKOrigin (size_t min[3]) const
 
void GetDimensions (size_t dims[3]) const
 
float GetMissingValue () const
 
void SetMissingValue (float missing_value)
 
bool HasMissingData () const
 
virtual int GetInterpolationOrder () const
 
virtual void SetInterpolationOrder (int order)
 
virtual void GetRange (float range[2]) const
 
virtual void HasPeriodic (bool *idim, bool *jdim, bool *kdim) const
 
void GetBlockSize (size_t bs[3]) const
 
virtual int GetRank ()
 
float ** GetBlks () const
 
size_t GetNumBlks () const
 
float Next ()
 
void ResetItr ()
 
Iterator begin ()
 
Iterator end ()
 
ConstIterator begin () const
 
ConstIterator end () const
 

Additional Inherited Members

- Protected Member Functions inherited from VAPoR::RegularGrid
float & _AccessIJK (float **blks, size_t i, size_t j, size_t k) const
 
void _ClampCoord (double &x, double &y, double &z) const
 
void _SetExtents (const double extents[6])
 

Detailed Description

Definition at line 25 of file LayeredGrid.h.

Constructor & Destructor Documentation

VAPoR::LayeredGrid::LayeredGrid ( const size_t  bs[3],
const size_t  min[3],
const size_t  max[3],
const double  extents[6],
const bool  periodic[3],
float **  blks,
float **  coords,
int  varying_dim 
)

Construct a layered grid sampling a 3D or 2D scalar function

Parameters
[in]bsA three-element vector specifying the dimensions of each block storing the sampled scalar function.
[in]minA three-element vector specifying the ijk index of the first point in the grid. The first grid point need not coincide with block boundaries. I.e. the indecies need not be (0,0,0)
[in]maxA three-element vector specifying the ijk index of the last point in the grid
[in]extentsA six-element vector specifying the user coordinates of the first (first three elements) and last (last three elements) of the grid points indicated by min and max, respectively. The elements of the extents parameter corresponding to the varying dimension are ignored. The varying dimension extents are provided by the coords parameter.
[in]periodicA three-element boolean vector indicating which i,j,k indecies, respectively, are periodic. The varying dimension may not be periodic.
[in]blksAn array of blocks containing the sampled function. The dimensions of each block is given by bs. The number of blocks is given by the product of the terms:
(max[i]/bs[i] - min[i]/bs[i] + 1)

over i = 0..2.

Parameters
[in]coordsAn array of blocks with dimension and number the same as blks specifying the varying dimension grid point coordinates.
[in]varying_dimAn enumerant indicating which axis is the varying dimension: 0 for I, 1 for J, 2 for K
VAPoR::LayeredGrid::LayeredGrid ( const size_t  bs[3],
const size_t  min[3],
const size_t  max[3],
const double  extents[6],
const bool  periodic[3],
float **  blks,
float **  coords,
int  varying_dim,
float  missing_value 
)

Construct a layered grid sampling a 3D or 2D scalar function that contains missing values.

This constructor adds a parameter, missing_value, that specifies the value of missing values in the sampled function. When reconstructing the function at arbitrary coordinates special consideration is given to grid points with missing values that are used in the reconstruction.

See also
GetValue()
virtual VAPoR::LayeredGrid::~LayeredGrid ( )
virtual

Member Function Documentation

virtual void VAPoR::LayeredGrid::GetBoundingBox ( const size_t  min[3],
const size_t  max[3],
double  extents[6] 
) const
virtual

Return the extents of the axis-aligned bounding box enclosign a region

This method returns min and max extents, in user coordinates, of the smallest axis-aligned box enclosing the region defined by the corner grid points, min and max.

Note
The results are undefined if any coordinate of min is greater than the coresponding coordinate of max.
Parameters
[out]extentsA six-element array, the first three values will contain the minimum coordinate, and the last three values the maximum coordinate
See also
GetDimensions(), RegularGrid()

Reimplemented from VAPoR::RegularGrid.

float** VAPoR::LayeredGrid::GetCoordBlks ( ) const
inline

Return the internal data structure containing a copy of the coordinate blocks passed in by the constructor

Definition at line 180 of file LayeredGrid.h.

virtual void VAPoR::LayeredGrid::GetEnclosingRegion ( const double  minu[3],
const double  maxu[3],
size_t  min[3],
size_t  max[3] 
) const
virtual

Get voxel coordinates of grid containing a region

Calculates the starting and ending IJK voxel coordinates of the smallest grid completely containing the rectangular region defined by the user coordinates minu and maxu If rectangluar region defined by minu and maxu can not be contained the minimum and maximum IJK coordinates are returned in min and max, respectively

Parameters
[in]minuUser coordinates of minimum coorner
[in]maxuUser coordinates of maximum coorner
[out]minInteger coordinates of minimum coorner
[out]maxInteger coordinates of maximum coorner

Reimplemented from VAPoR::RegularGrid.

void VAPoR::LayeredGrid::GetIJKIndex ( double  x,
double  y,
double  z,
size_t *  i,
size_t *  j,
size_t *  k 
) const
virtual

Return the closest grid point to the specified user coordinates

This method returns the ijk index of the grid point closest to the specified user coordinates based on Euclidean distance. If any of the input coordinates correspond to periodic dimensions the the coordinate(s) are first re-mapped to lie inside the grid extents as returned by GetUserExtents()

Parameters
[in]xcoordinate along fastest varying dimension
[in]ycoordinate along second fastest varying dimension
[in]zcoordinate along third fastest varying dimension
[out]iindex of grid point along fastest varying dimension
[out]jindex of grid point along second fastest varying dimension
[out]kindex of grid point along third fastest varying dimension

Reimplemented from VAPoR::RegularGrid.

void VAPoR::LayeredGrid::GetIJKIndexFloor ( double  x,
double  y,
double  z,
size_t *  i,
size_t *  j,
size_t *  k 
) const
virtual

Return the corner grid point of the cell containing the specified user coordinates

This method returns the smallest ijk index of the grid point of associated with the cell containing the specified user coordinates. If any of the input coordinates correspond to periodic dimensions the the coordinate(s) are first re-mapped to lie inside the grid extents as returned by GetUserExtents()

If the specified coordinates lie outside of the grid (are not contained by any cell) the lowest valued ijk index of the grid points defining the boundary cell closest to the point are returned.

Parameters
[in]xcoordinate along fastest varying dimension
[in]ycoordinate along second fastest varying dimension
[in]zcoordinate along third fastest varying dimension
[out]iindex of grid point along fastest varying dimension
[out]jindex of grid point along second fastest varying dimension
[out]kindex of grid point along third fastest varying dimension

Reimplemented from VAPoR::RegularGrid.

virtual void VAPoR::LayeredGrid::GetMinCellExtents ( double *  x,
double *  y,
double *  z 
) const
virtual

Return the minimum grid spacing between all grid points

This method returns the minimum distance, in user coordinates, between adjacent grid points for all cells in the grid.

Parameters
[out]xMinimum distance between grid points along X axis
[out]yMinimum distance between grid points along Y axis
[out]zMinimum distance between grid points along Z axis
Note
For a regular grid all cells have the same dimensions

Reimplemented from VAPoR::RegularGrid.

int VAPoR::LayeredGrid::GetUserCoordinates ( size_t  i,
size_t  j,
size_t  k,
double *  x,
double *  y,
double *  z 
) const
virtual

Return the user coordinates of a grid point

This method returns the user coordinates of the grid point specified by index(i,j,k)

Parameters
[in]iindex of grid point along fastest varying dimension
[in]jindex of grid point along second fastest varying dimension
[in]kindex of grid point along third fastest varying dimension
[out]xcoordinate of grid point along fastest varying dimension
[out]ycoordinate of grid point along second fastest varying dimension
[out]zcoordinate of grid point along third fastest varying dimension
Return values
statusA negative int is returned if index(i,j,k) is out out of bounds

Reimplemented from VAPoR::RegularGrid.

virtual void VAPoR::LayeredGrid::GetUserExtents ( double  extents[6]) const
inlinevirtual

Return the extents of the user coordinate system

This method returns min and max extents of the user coordinate system defined on the grid. The minimum extent is the coordinate of the grid point with index(0,0,0). The maximum extent is the coordinate of the grid point with index(nx-1, ny-1, nz-1), where nx, ny, and nz are the I, J, K dimensions, respectively.

Parameters
[out]extentsA six-element array, the first three values will contain the minimum coordinate, and the last three values the maximum coordinate
See also
GetDimensions(), RegularGrid()

Reimplemented from VAPoR::RegularGrid.

Definition at line 106 of file LayeredGrid.h.

float VAPoR::LayeredGrid::GetValue ( double  x,
double  y,
double  z 
) const
virtual

Get the reconstructed value of the sampled scalar function

This method reconstructs the scalar field at an arbitrary point in space. If the point's coordinates are outside of the grid's coordinate extents as returned by GetUserExtents(), and the grid is not periodic along the out-of-bounds axis, the value returned will be the missing_value.

If the value of any of the grid point samples used in the reconstruction is the missing_value then the result returned is the missing_value.

The reconstruction method used is determined by interpolation order returned by GetInterpolationOrder()

Parameters
[in]xcoordinate along fastest varying dimension
[in]ycoordinate along second fastest varying dimension
[in]zcoordinate along third fastest varying dimension
See also
GetInterpolationOrder(), HasPeriodic(), GetMissingValue()
GetUserExtents()

Reimplemented from VAPoR::RegularGrid.

bool VAPoR::LayeredGrid::InsideGrid ( double  x,
double  y,
double  z 
) const
virtual

Return true if the specified point lies inside the grid

This method can be used to determine if a point expressed in user coordinates reside inside or outside the grid

Parameters
[in]xcoordinate along fastest varying dimension
[in]ycoordinate along second fastest varying dimension
[in]zcoordinate along third fastest varying dimension
Return values
boolTrue if point is inside the grid

Reimplemented from VAPoR::RegularGrid.

int VAPoR::LayeredGrid::Reshape ( const size_t  min[3],
const size_t  max[3],
const bool  periodic[3] 
)
virtual

Change the voxel exents specified by the constructor

This method permits the grid to be reshaped under the constraints that: 1) the number of blocks does not change, and 2) the grid can only get smaller, not larger. I.e. min can only increase and max can only decrease.

Parameters
[in]minA three-element vector specifying the ijk index of the first point in the grid. The first grid point need not coincide with block boundaries. I.e. the indecies need not be (0,0,0): the first grid point is not required to be the first element of the array
[in]maxA three-element vector specifying the ijk index of the last point in the grid.
[in]periodicA three-element boolean vector indicating which i,j,k indecies, respectively, are periodic

Reimplemented from VAPoR::RegularGrid.

virtual void VAPoR::LayeredGrid::SetPeriodic ( const bool  periodic[3])
virtual

Set periodic boundaries

This method changes the periodicity of boundaries set by the class constructor. It overides the base class method to ensure the varying dimension is not periodic.

Parameters
[in]periodicA three-element boolean vector indicating which i,j,k indecies, respectively, are periodic. The varying dimension may not be periodic.

Reimplemented from VAPoR::RegularGrid.


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