VAPOR3 3.9.4
|
Abstract base class for a 2D or 3D structured or unstructured grid. More...
#include <Grid.h>
Classes | |
class | AbstractIterator |
class | ConstCellIteratorBoxSG |
class | ConstCellIteratorSG |
class | ConstNodeIteratorBoxSG |
class | ConstNodeIteratorSG |
class | ForwardIterator |
class | InsideBox |
class | PolyIterator |
Public Types | |
typedef const CoordType | ConstCoordType |
typedef Grid::PolyIterator< ConstCoordType > | ConstCoordItr |
typedef Grid::AbstractIterator< ConstCoordType > | ConstCoordItrAbstract |
typedef const DimsType | ConstIndexType |
typedef Grid::PolyIterator< ConstIndexType > | ConstNodeIterator |
typedef Grid::AbstractIterator< ConstIndexType > | ConstNodeIteratorAbstract |
typedef Grid::PolyIterator< ConstIndexType > | ConstCellIterator |
Cell index iterator. Iterates over grid cell indices. | |
typedef Grid::AbstractIterator< ConstIndexType > | ConstCellIteratorAbstract |
typedef Grid::ForwardIterator< Grid > | Iterator |
typedef Grid::ForwardIterator< Grid const > | ConstIterator |
Public Member Functions | |
Grid (const std::vector< size_t > &dims, const std::vector< size_t > &bs, const std::vector< float * > &blks, size_t topology_dimension) | |
Grid (const DimsType &dims, const DimsType &bs, const std::vector< float * > &blks, size_t topology_dimension) | |
Grid () | |
virtual | ~Grid ()=default |
const DimsType & | GetDimensions () const |
size_t | GetNumDimensions () const |
virtual DimsType | GetCoordDimensions (size_t dim) const =0 |
virtual std::string | GetType () const =0 |
virtual size_t | GetGeometryDim () const =0 |
virtual const DimsType & | GetNodeDimensions () const =0 |
virtual const size_t | GetNumNodeDimensions () const =0 |
virtual const DimsType & | GetCellDimensions () const =0 |
virtual const size_t | GetNumCellDimensions () const =0 |
const std::vector< size_t > | GetDimensionInBlks () const |
const std::vector< size_t > & | GetBlockSize () const |
const std::vector< float * > & | GetBlks () const |
virtual float | GetValueAtIndex (const DimsType &indices) const |
virtual float | GetValueAtIndex (const std::vector< size_t > &indices) const |
virtual void | SetValue (const DimsType &indices, float value) |
virtual void | SetValue (const size_t indices[3], float value) |
virtual float | AccessIJK (size_t i, size_t j=0, size_t k=0) const |
void | SetValueIJK (size_t i, size_t j, size_t k, float v) |
void | SetValueIJK (size_t i, size_t j, float v) |
void | SetValueIJK (size_t i, float v) |
virtual float | GetValue (const CoordType &coords) const |
virtual float | GetValue (const std::vector< double > &coords) const |
virtual float | GetValue (const double coords[]) const |
virtual float | GetValue (double x, double y) const |
virtual float | GetValue (double x, double y, double z) const |
virtual void | GetUserExtents (CoordType &minu, CoordType &maxu) const |
virtual void | GetUserExtents (double minu[3], double maxu[3]) const |
virtual void | GetUserExtents (std::vector< double > &minu, std::vector< double > &maxu) const |
virtual void | GetBoundingBox (const DimsType &min, const DimsType &max, CoordType &minu, CoordType &maxu) const =0 |
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 bool | GetEnclosingRegion (const CoordType &minu, const CoordType &maxu, DimsType &min, DimsType &max) const =0 |
virtual float | GetMissingValue () const |
void | SetMissingValue (float missing_value) |
void | SetHasMissingValues (bool flag) |
bool | HasMissingData () const |
virtual bool | HasInvertedCoordinateSystemHandiness () const |
virtual int | GetInterpolationOrder () const |
virtual void | SetInterpolationOrder (int order) |
virtual void | GetUserCoordinates (const DimsType &indices, CoordType &coords) const =0 |
virtual void | GetUserCoordinates (const size_t indices[], double coords[]) const |
virtual void | GetUserCoordinates (const std::vector< size_t > &indices, std::vector< double > &coords) const |
virtual void | GetUserCoordinates (size_t i, double &x, double &y, double &z) const |
virtual void | GetUserCoordinates (size_t i, size_t j, double &x, double &y, double &z) const |
virtual void | GetUserCoordinates (size_t i, size_t j, size_t k, double &x, double &y, double &z) const |
virtual bool | GetIndicesCell (const CoordType &coords, DimsType &indices) const =0 |
virtual bool | GetIndicesCell (const double coords[3], size_t indices[3]) const |
virtual bool | GetIndicesCell (const std::vector< double > &coords, std::vector< size_t > &indices) const |
virtual void | GetRange (float range[2]) const |
virtual void | GetRange (const DimsType &min, const DimsType &max, float range[2]) const |
virtual void | GetRange (std::vector< size_t > min, std::vector< size_t > max, float range[2]) const |
virtual bool | InsideGrid (const CoordType &coords) const =0 |
virtual bool | InsideGrid (const double coords[3]) const |
virtual bool | InsideGrid (const std::vector< double > &coords) const |
virtual bool | GetCellNodes (const DimsType &cindices, std::vector< DimsType > &nodes) const =0 |
virtual bool | GetCellNodes (const size_t cindices[], std::vector< DimsType > &nodes) const |
virtual bool | GetCellNeighbors (const DimsType &cindices, std::vector< DimsType > &cells) const =0 |
virtual bool | GetNodeCells (const DimsType &indices, std::vector< DimsType > &cells) const =0 |
virtual size_t | GetMaxVertexPerFace () const =0 |
virtual size_t | GetMaxVertexPerCell () const =0 |
virtual void | ClampCoord (const CoordType &coords, CoordType &cCoords) const =0 |
virtual void | ClampCoord (const double coords[3], double cCoords[3]) const |
virtual void | ClampIndex (const DimsType &indices, DimsType &cIndices) const |
virtual void | ClampCellIndex (const DimsType &indices, DimsType &cIndices) const |
virtual void | SetPeriodic (const std::vector< bool > &periodic) |
virtual const std::vector< bool > & | GetPeriodic () const |
virtual size_t | GetTopologyDim () const |
virtual long | GetNodeOffset () const |
virtual void | SetNodeOffset (long offset) |
virtual long | GetCellOffset () const |
virtual void | SetCellOffset (long offset) |
virtual DimsType | GetMinAbs () const |
virtual void | SetMinAbs (const DimsType &minAbs) |
virtual ConstCoordItr | ConstCoordBegin () const =0 |
Return constant grid coordinate iterator. | |
virtual ConstCoordItr | ConstCoordEnd () const =0 |
virtual ConstNodeIterator | ConstNodeBegin () const |
virtual ConstNodeIterator | ConstNodeBegin (const CoordType &minu, const CoordType &maxu) const |
virtual ConstNodeIterator | ConstNodeBegin (const std::vector< double > &minu, const std::vector< double > &maxu) const |
virtual ConstNodeIterator | ConstNodeEnd () const |
virtual ConstCellIterator | ConstCellBegin () const |
virtual ConstCellIterator | ConstCellBegin (const CoordType &minu, const CoordType &maxu) const |
virtual ConstCellIterator | ConstCellBegin (const std::vector< double > &minu, const std::vector< double > &maxu) const |
virtual ConstCellIterator | ConstCellEnd () const |
Iterator | begin (const CoordType &minu, const CoordType &maxu) |
Iterator | begin (const std::vector< double > &minu, const std::vector< double > &maxu) |
Iterator | begin () |
Iterator | end () |
ConstIterator | cbegin (const CoordType &minu, const CoordType &maxu) const |
ConstIterator | cbegin (const std::vector< double > &minu, const std::vector< double > &maxu) |
ConstIterator | cbegin () const |
ConstIterator | cend () const |
Static Public Member Functions | |
static size_t | GetNumDimensions (DimsType dims) |
static DimsType | Dims (const DimsType &min, const DimsType &max) |
static bool | PointInsideBoundingRectangle (const double pt[], const double verts[], int n) |
template<typename T > | |
static void | CopyToArr3 (const std::vector< T > &src, std::array< T, 3 > &dst) |
template<typename T > | |
static void | CopyToArr3 (const T *src, size_t n, std::array< T, 3 > &dst) |
template<typename T > | |
static void | CopyFromArr3 (const std::array< T, 3 > &src, std::vector< T > &dst) |
template<typename T > | |
static void | CopyFromArr3 (const std::array< T, 3 > &src, T *dst) |
Protected Member Functions | |
virtual float | GetValueNearestNeighbor (const CoordType &coords) const =0 |
virtual float | GetValueLinear (const CoordType &coords) const =0 |
virtual void | GetUserExtentsHelper (CoordType &minu, CoordType &maxu) const =0 |
virtual float * | GetValuePtrAtIndex (const std::vector< float * > &blks, const DimsType &indices) const |
virtual void | ClampIndex (const std::vector< size_t > &dims, const DimsType indices, DimsType &cIndices) const |
virtual void | ClampIndex (const DimsType &dims, const DimsType indices, DimsType &cIndices) const |
float | BilinearInterpolate (size_t i, size_t j, size_t k, const double xwgt, const double ywgt) const |
float | TrilinearInterpolate (size_t i, size_t j, size_t k, const double xwgt, const double ywgt, const double zwgt) const |
Friends | |
VDF_API friend std::ostream & | operator<< (std::ostream &o, const Grid &g) |
Abstract base class for a 2D or 3D structured or unstructured grid.
This abstract base class defines a 2D or 3D structured or unstructured grid.
The grid samples a scalar function at each grid point. Each grid point can be addressed by a 3-element array defined as type DimsType
. The fastest varying dimension is given by indices[0], etc.
Because grid samples are repesented internally as arrays, when accessing multiple grid points better performance is achieved by using unit stride.
indices | A DimsType i in the range 0..max, where max is one minus the value of the corresponding element returned by GetDimensions(). |
coords | A CoordType containig the coordinates of a point in user-defined coordinates. Elements with indices greater than GetGeometryDim()-1 are ignored. |
typedef const CoordType VAPoR::Grid::ConstCoordType |
Coordinate iterator. Iterates over grid node/cell coordinates
The ConstCoordItr can be dererenced to return a grid node or cell's coordinates. The determination, node or cell, is determined by the location of the sampled data within the grid (node, face, cell, etc)
N.B. Current only works with node coordinates
typedef const DimsType VAPoR::Grid::ConstIndexType |
typedef Grid::ForwardIterator<Grid const> VAPoR::Grid::ConstIterator |
VAPoR::Grid::Grid | ( | const std::vector< size_t > & | dims, |
const std::vector< size_t > & | bs, | ||
const std::vector< float * > & | blks, | ||
size_t | topology_dimension | ||
) |
Construct a structured or unstructured grid sampling a 3D or 2D scalar function
[in] | dims | Dimensions of arrays containing grid data. |
[in] | bs | A DimsType with, specifying the dimensions of each block storing the sampled scalar function. |
[in] | blks | An 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: |
over i = 0..dims.size()-1
A shallow copy of the blocks is made by the constructor. Memory referenced by the elements of blks
should remain valid until the class instance is destroyed.
[in] | topology_dimension | Topological dimension of grid: 1, 2 or 3, for 1D, 2D or 3D, respectively. Grids with 2D topology are composed of 2D polygons, while grids with 3D topology are composed of 3D polyhedron |
VAPoR::Grid::Grid | ( | const DimsType & | dims, |
const DimsType & | bs, | ||
const std::vector< float * > & | blks, | ||
size_t | topology_dimension | ||
) |
VAPoR::Grid::Grid | ( | ) |
|
virtualdefault |
|
virtual |
This method provides an alternate interface to Grid::GetValueAtIndex() If the dimensionality of the grid as determined by GetDimensions() is less than three subsequent parameters are ignored. Parameters that are outside of range are clamped to boundaries.
[in] | i | Index into first fastest varying dimension |
[in] | j | Index into second fastest varying dimension |
[in] | k | Index into third fastest varying dimension |
|
inline |
|
protected |
|
inline |
|
inline |
|
inline |
|
inline |
|
inlinevirtual |
Clamp grid cell indices
Same as ClampIndex() accept that indices are clamped to GetCellDimensions()
|
pure virtual |
Clamp periodic coordinates and ensure valid coordinate vector dimension
This method ensures that periodic coordinates are within the bounding box of the grid and that the coordinate vector dimension does not exceed the number of allowable coordinates as returned by GetGeometryDim().
[in] | coords | A coordinate vector |
[out] | cCoords | The clamped coordintes coords |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
inlinevirtual |
Reimplemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
inlinevirtual |
Clamp grid array indices
This method ensures that grid indices are not out of bounds, clamping any elements of indices
that exceeds or equals the corresponding element of GetNodeDimensions() to that value minus 1.
[in] | indices | An array index vector |
[out] | cindices | An array index vector, clamped as described above |
Definition at line 688 of file Grid.h.
References ClampIndex().
Referenced by ClampIndex().
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
pure virtual |
Return constant grid coordinate iterator.
Implemented in VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
pure virtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
Compute the dimensions of a rectangular region bounded by min
and max
coordinates
|
inline |
|
inline |
|
pure virtual |
Return the extents of the axis-aligned bounding box enclosign a region
This pure virtual method returns min and max extents, in user coordinates, of the smallest axis-aligned box enclosing the region defined by the grid indices, min
and max
. Every grid point in the range min
to max
will be contained in, or reside on, the box (rectangle) whose extents are given by minu
and maxu
.
min
is greater than the coresponding coordinate of max
, or if max
is outside of the valid dimension range (See GetDimension()).The size of min
and max
must match the grid's dimension as returned by GetDimension()
[in] | min | An array containing the minimum grid indices (offsets). |
[in] | max | An array containing the maximum grid indices (offsets). |
[out] | minu | A two-to-three element array containing the minimum user coordinates. |
[out] | maxu | A two-to-three element array containing the maximum user coordinates. |
Implemented in VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
inlinevirtual |
|
pure virtual |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
pure virtual |
Get the IDs (indices) of all of the cells that border a cell
This method returns a vector of index vectors. Each index vector contains the indices of a cell that borders the cell given by indices
. If a cell edge is a boundary edge, having no neighbors, the associated index vector for that border will be empty. The cell IDs are returned in counter-clockwise order
[in] | cindices | An ordered vector of multi-dimensional cell indices. |
[out] | cells | A vector of index vectors. Each index vector has size given by GetDimensions.size() |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
pure virtual |
Get the indices of the nodes that define a cell
This method returns a vector of index vectors. Each index vector contains the indices for a node that defines the given cell ID indices
.
For 2D cells the node IDs are returned in counter-clockwise order. For 3D cells the ordering is dependent on the shape of the node. TBD.
[in] | cindices | An array of size Grid::GetDimensions.size() specifying the index of the cell to query. |
[out] | nodes | An array of size Grid::GetMaxVertexPerCell() x Grid::GetDimensions.size() that will contain the concatenated list of node indices on return. |
[out] | n | The number of node indices returned in nodes , an integer in the range (0..Grid::GetMaxVertexPerCell()). |
Implemented in VAPoR::StructuredGrid, VAPoR::UnstructuredGrid, and VAPoR::StructuredGrid.
|
inlinevirtual |
|
inlinevirtual |
|
pure virtual |
Return the dimensions of the specified coordinate variable
[in] | dim | An integer between 0 and the return value of GetGeometryDim() - 1, indicating which coordinate dimensions are to be returned. |
if dim
is greater than or equal to GetGeometryDim() - 1 a vector of length one with its single component equal to one is returned.
Implemented in VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
inline |
Return the ijk dimensions of grid in blocks
Returns the number of blocks defined along each axis of the grid
[out] | dims | A two or three element array containing the grid dimension in blocks |
|
inline |
Return the dimensions of grid connectivity array
[out] | dims | The value of dims parameter provided to the constructor. If the parameter has less than 3 values, then number 1 will be filled. |
|
pure virtual |
Get bounding indices of grid containing a region
Calculates the starting and ending grid indices 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 indices are returned in min
and max
, respectively
[in] | minu | User coordinates of minimum coorner |
[in] | maxu | User coordinates of maximum coorner |
[out] | min | Integer coordinates of minimum coorner |
[out] | max | Integer coordinates of maximum coorner |
Implemented in VAPoR::StructuredGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
pure virtual |
Get geometric dimension of cells
Returns the geometric dimension of the cells in the mesh. I.e. the number of spatial coordinates for each grid point. Valid values are 0..3. The geometric dimension must be equal to or greater than the topology dimension.
Implemented in VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
pure virtual |
Return the indices of the cell containing the specified user coordinates
This method returns the cell ID (index) of 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 method returns false, otherwise true is returned.
status | true on success, false if the point is not contained by any cell. |
Implemented in VAPoR::UnstructuredGridCoordless, VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridLayered, VAPoR::LayeredGrid, and VAPoR::UnstructuredGrid2D.
|
inlinevirtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::UnstructuredGrid2D.
|
inlinevirtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::UnstructuredGrid2D.
|
inlinevirtual |
Return the interpolation order to be used during function reconstruction
This method returns the order of the interpolation method that will be used when reconstructing the sampled scalar function
Reimplemented in VAPoR::LayeredGrid.
|
pure virtual |
Return the maximum number of vertices per cell
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
pure virtual |
Return the maximum number of vertices per cell face
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
inlinevirtual |
Return the absolute minimum grid coordinate
This method returns the absolute minimum grid coordinate. If the Grid contains a subregion extracted from a larger mesh this absolute minimum grid coordinate gives the offset of the first gridpoint in this grid relative to the larger mesh.
|
virtual |
Return the value of the missing_value parameter
The missing value is a special value intended to indicate that the value of the sampled or reconstructed function is unknown at a particular point.
|
pure virtual |
Get the IDs (indices) of all of the cells that share a node
This method returns a vector of index vectors. Each index vector contains the indices of a cell that share the node given by indices
. The cell IDs are returned in counter-clockwise order
[out] | nodes | A vector of index vectors . Each index vector has size given by GetDimensions.size() |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
pure virtual |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
inlinevirtual |
|
pure virtual |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
inline |
Return the useful number of dimensions of grid connectivity array
[out] | dims | The number of non-unit components of dims parameter provided to the constructor. |
Definition at line 107 of file Grid.h.
References GetNumDimensions().
Referenced by GetNumDimensions().
|
static |
Return the number of non-unit dimensions
Return the number of non-unit dimensions in dims
[in] | dims |
|
pure virtual |
Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.
|
inlinevirtual |
|
virtual |
|
virtual |
Return the min and max data value
This method returns the values of grid points with min and max values, respectively.
For dataless grids the missing value is returned.
[out] | range[2] | A two-element array containing the mininum and maximum values, in that order |
|
inlinevirtual |
|
inlinevirtual |
Get topological dimension of the mesh
Return the number of topological dimensions for the mesh cells. Valid values are in the range 0..3, with, for example, 0 corresponding to points; 1 to lines; 2 to triangles, quadrilaterals, etc.; and 3 to hexahedron, tetrahedron, etc.
Reimplemented in VAPoR::ConstantGrid.
|
pure virtual |
Implemented in VAPoR::ArbitrarilyOrientedRegularGrid, VAPoR::ConstantGrid, VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::StructuredGrid, VAPoR::UnstructuredGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
pure virtual |
Return the user coordinates of a grid point
This method returns the user coordinates of the grid point specified by indices
Results are undefined if indices
is out of range.
[in] | index | Logical index into coordinate array. The dimensionality and range of component values are given by GetNodeDimensions(). The starting value for each component of index is zero. index must contain GetNodeDimensions().size() elements. |
[out] | coord | User coordinates of grid point with indices given by indices . coord must have space for the number of elements returned by GetGeometryDim(). |
Implemented in VAPoR::ArbitrarilyOrientedRegularGrid, VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, VAPoR::UnstructuredGridLayered, VAPoR::LayeredGrid, and VAPoR::RegularGrid.
|
inlinevirtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.
|
inlinevirtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.
|
virtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.
|
virtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.
|
virtual |
Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.
Return the extents of the user coordinate system
This pure virtual method returns min and max extents of the user coordinate system defined on the grid. The extents of the returned box are guaranteed to contain all points in the grid.
[out] | minu | A two or three element array containing the minimum user coordinates. |
[out] | maxu | A two or three element array containing the maximum user coordinates. |
|
inlinevirtual |
|
inlinevirtual |
|
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()
[in] | coords | A vector of size matching the topology dimension of the mesh whose contents specify the coordinates of a point in space. |
Reimplemented in VAPoR::ConstantGrid, and VAPoR::LayeredGrid.
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
Reimplemented in VAPoR::SphericalGrid.
|
virtual |
Get the data value at the indicated grid point
This method provides read access to the scalar data value defined at the grid point indicated by indices
. The range of valid indices is between zero and dim - 1, where dim is the dimesion of the grid returned by GetDimensions()
If any of the indecies
are outside of the valid range the results are undefined
[in] | indices | of grid point along fastest varying dimension. |
|
inlinevirtual |
|
protectedpure virtual |
|
protectedpure virtual |
|
protectedvirtual |
|
inlinevirtual |
Return true if mesh primitives have counter clockwise winding order.
Reimplemented in VAPoR::StructuredGrid.
|
inline |
Return missing data flag
This method returns true if the missing data flag is set. This does not imply that grid points exist with missing data, only that the the grid points should be compared against the value of GetMissingValue() whenever operations are performed on them.
|
pure 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
[in] | coords | User coordinates of point in space |
bool | True if point is inside the grid |
Implemented in VAPoR::ConstantGrid, VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.
|
inlinevirtual |
|
inlinevirtual |
|
inlinestatic |
Test whether a point is contained in a bounding rectangle
This static method checks to see if a 2D point pt
is contained in the smallest rectangle that bounds the list of 2D points (vertices) given by verts
. If pt
is inside or on the boundary of the bounding rectangle true is returned, otherwise false
[in] | pt | A two-element array of point coordinates |
[in] | verts | An array of dimensions n * 2 containing a list of points . |
[in] | n | The number of 2D points in verts |
status | Returns true if \pt is inside or on the bounding rectangle of verts . False otherwise. |
Definition at line 781 of file Grid.h.
References VAssert.
|
inlinevirtual |
Reimplemented in VAPoR::UnstructuredGrid.
|
inline |
Set missing data flag
Set this flag if missing values may exist in the data. This missing values are denoted by the value GetMissingValue(). Subsequent processing of grid data will be compared against the value of GetMissingValue() if this flag is set.
|
virtual |
Set the interpolation order to be used during function reconstruction
This method sets the order of the interpolation method that will be used when reconstructing the sampled scalar function. Valid values of order
are 0 and 1, corresponding to nearest-neighbor and linear interpolation, respectively. If order
is invalid it will be silently set to 1. The default interpolation order is 1
[in] | order | interpolation order |
Reimplemented in VAPoR::LayeredGrid.
|
inlinevirtual |
Set the absolute minimum grid coordinate
[in] | minAbs | Must have same dimensionality as constructors dims parameter. Otherwise may contain any value, but is intended to contain the offset to the first grid point in the mesh. The default is the zero vector |
|
inline |
Set the missing value indicator
This method sets the value of the missing value indicator. The method does not alter the value of any grid point locations, nor does it alter the missing data flag set with the constructor.
|
inlinevirtual |
Reimplemented in VAPoR::UnstructuredGrid.
|
inlinevirtual |
Set periodic boundaries
This method changes the periodicity of boundaries set by the class constructor
[in] | periodic | A boolean vector of size given by GetGeometryDim() indicating which coordinate axes are periodic. |
Reimplemented in VAPoR::LayeredGrid.
|
virtual |
Set the data value at the indicated grid point
This method sets the data value of the grid point indexed by indices
to value
.
|
inlinevirtual |
|
inline |
|
inline |
void VAPoR::Grid::SetValueIJK | ( | size_t | i, |
size_t | j, | ||
size_t | k, | ||
float | v | ||
) |
|
protected |