VAPOR3 3.9.4
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Friends | List of all members
VAPoR::Grid Class Referenceabstract

Abstract base class for a 2D or 3D structured or unstructured grid. More...

#include <Grid.h>

Inheritance diagram for VAPoR::Grid:
VAPoR::ConstantGrid VAPoR::StructuredGrid VAPoR::UnstructuredGrid VAPoR::CurvilinearGrid VAPoR::LayeredGrid VAPoR::RegularGrid VAPoR::StretchedGrid VAPoR::UnstructuredGrid2D VAPoR::UnstructuredGrid3D VAPoR::UnstructuredGridCoordless VAPoR::UnstructuredGridLayered

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< ConstCoordTypeConstCoordItr
 
typedef Grid::AbstractIterator< ConstCoordTypeConstCoordItrAbstract
 
typedef const DimsType ConstIndexType
 
typedef Grid::PolyIterator< ConstIndexTypeConstNodeIterator
 
typedef Grid::AbstractIterator< ConstIndexTypeConstNodeIteratorAbstract
 
typedef Grid::PolyIterator< ConstIndexTypeConstCellIterator
 Cell index iterator. Iterates over grid cell indices.
 
typedef Grid::AbstractIterator< ConstIndexTypeConstCellIteratorAbstract
 
typedef Grid::ForwardIterator< GridIterator
 
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 DimsTypeGetDimensions () 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 DimsTypeGetNodeDimensions () const =0
 
virtual const size_t GetNumNodeDimensions () const =0
 
virtual const DimsTypeGetCellDimensions () 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)
 

Detailed Description

Abstract base class for a 2D or 3D structured or unstructured grid.

Author
John Clyne

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.

Parameters
indicesA DimsType i in the range 0..max, where max is one minus the value of the corresponding element returned by GetDimensions().
coordsA CoordType containig the coordinates of a point in user-defined coordinates. Elements with indices greater than GetGeometryDim()-1 are ignored.

Definition at line 56 of file Grid.h.

Member Typedef Documentation

◆ ConstCellIterator

Cell index iterator. Iterates over grid cell indices.

Definition at line 956 of file Grid.h.

◆ ConstCellIteratorAbstract

Definition at line 957 of file Grid.h.

◆ ConstCoordItr

Definition at line 937 of file Grid.h.

◆ ConstCoordItrAbstract

Definition at line 938 of file Grid.h.

◆ 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

Definition at line 936 of file Grid.h.

◆ ConstIndexType

Node index iterator. Iterates over grid node indices

The ConstNodeIterator is dereferenced to give the index of a node within the grid

Definition at line 950 of file Grid.h.

◆ ConstIterator

Definition at line 1184 of file Grid.h.

◆ ConstNodeIterator

Definition at line 951 of file Grid.h.

◆ ConstNodeIteratorAbstract

Definition at line 952 of file Grid.h.

◆ Iterator

Definition at line 1183 of file Grid.h.

Constructor & Destructor Documentation

◆ Grid() [1/3]

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

Parameters
[in]dimsDimensions of arrays containing grid data.
[in]bsA DimsType with, specifying the dimensions of each block storing the sampled scalar function.
[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:
(((dims[i]-1) / bs[i]) + 1)

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.

Parameters
[in]topology_dimensionTopological 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

◆ Grid() [2/3]

VAPoR::Grid::Grid ( const DimsType dims,
const DimsType bs,
const std::vector< float * > &  blks,
size_t  topology_dimension 
)

◆ Grid() [3/3]

VAPoR::Grid::Grid ( )

◆ ~Grid()

virtual VAPoR::Grid::~Grid ( )
virtualdefault

Member Function Documentation

◆ AccessIJK()

virtual float VAPoR::Grid::AccessIJK ( size_t  i,
size_t  j = 0,
size_t  k = 0 
) const
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.

Parameters
[in]iIndex into first fastest varying dimension
[in]jIndex into second fastest varying dimension
[in]kIndex into third fastest varying dimension

◆ begin() [1/3]

Iterator VAPoR::Grid::begin ( )
inline

Definition at line 1199 of file Grid.h.

◆ begin() [2/3]

Iterator VAPoR::Grid::begin ( const CoordType minu,
const CoordType maxu 
)
inline

Construct a begin iterator that will iterate through elements inside or on the box defined by minu and maxu

Definition at line 1189 of file Grid.h.

◆ begin() [3/3]

Iterator VAPoR::Grid::begin ( const std::vector< double > &  minu,
const std::vector< double > &  maxu 
)
inline

Definition at line 1191 of file Grid.h.

◆ BilinearInterpolate()

float VAPoR::Grid::BilinearInterpolate ( size_t  i,
size_t  j,
size_t  k,
const double  xwgt,
const double  ywgt 
) const
protected

◆ cbegin() [1/3]

ConstIterator VAPoR::Grid::cbegin ( ) const
inline

Definition at line 1212 of file Grid.h.

◆ cbegin() [2/3]

ConstIterator VAPoR::Grid::cbegin ( const CoordType minu,
const CoordType maxu 
) const
inline

Definition at line 1203 of file Grid.h.

◆ cbegin() [3/3]

ConstIterator VAPoR::Grid::cbegin ( const std::vector< double > &  minu,
const std::vector< double > &  maxu 
)
inline

Definition at line 1205 of file Grid.h.

◆ cend()

ConstIterator VAPoR::Grid::cend ( ) const
inline

Definition at line 1214 of file Grid.h.

◆ ClampCellIndex()

virtual void VAPoR::Grid::ClampCellIndex ( const DimsType indices,
DimsType cIndices 
) const
inlinevirtual

Clamp grid cell indices

Same as ClampIndex() accept that indices are clamped to GetCellDimensions()

See also
ClampIndex()

Definition at line 696 of file Grid.h.

◆ ClampCoord() [1/2]

virtual void VAPoR::Grid::ClampCoord ( const CoordType coords,
CoordType cCoords 
) const
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().

Parameters
[in]coordsA coordinate vector
[out]cCoordsThe clamped coordintes coords
See also
GetGeometryDim()

Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.

◆ ClampCoord() [2/2]

virtual void VAPoR::Grid::ClampCoord ( const double  coords[3],
double  cCoords[3] 
) const
inlinevirtual
Deprecated:

Reimplemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.

Definition at line 667 of file Grid.h.

◆ ClampIndex() [1/3]

virtual void VAPoR::Grid::ClampIndex ( const DimsType dims,
const DimsType  indices,
DimsType cIndices 
) const
inlineprotectedvirtual

Definition at line 1255 of file Grid.h.

◆ ClampIndex() [2/3]

virtual void VAPoR::Grid::ClampIndex ( const DimsType indices,
DimsType cIndices 
) const
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.

Parameters
[in]indicesAn array index vector
[out]cindicesAn array index vector, clamped as described above
See also
GetNodeDimensions()

Definition at line 688 of file Grid.h.

References ClampIndex().

Referenced by ClampIndex().

◆ ClampIndex() [3/3]

virtual void VAPoR::Grid::ClampIndex ( const std::vector< size_t > &  dims,
const DimsType  indices,
DimsType cIndices 
) const
inlineprotectedvirtual

Definition at line 1245 of file Grid.h.

◆ ConstCellBegin() [1/3]

virtual ConstCellIterator VAPoR::Grid::ConstCellBegin ( ) const
inlinevirtual

Return constant grid cell coordinate iterator

If minu and maxu are specified the iterator is constrained to operation within the axis-aligned box defined by minu and maxu.

Parameters
[in]minuMinimum box coordinate.
[in]maxuMaximum box coordinate.

Definition at line 1096 of file Grid.h.

◆ ConstCellBegin() [2/3]

virtual ConstCellIterator VAPoR::Grid::ConstCellBegin ( const CoordType minu,
const CoordType maxu 
) const
inlinevirtual

Definition at line 1098 of file Grid.h.

◆ ConstCellBegin() [3/3]

virtual ConstCellIterator VAPoR::Grid::ConstCellBegin ( const std::vector< double > &  minu,
const std::vector< double > &  maxu 
) const
inlinevirtual

Definition at line 1103 of file Grid.h.

◆ ConstCellEnd()

virtual ConstCellIterator VAPoR::Grid::ConstCellEnd ( ) const
inlinevirtual

Definition at line 1111 of file Grid.h.

◆ ConstCoordBegin()

virtual ConstCoordItr VAPoR::Grid::ConstCoordBegin ( ) const
pure virtual

◆ ConstCoordEnd()

virtual ConstCoordItr VAPoR::Grid::ConstCoordEnd ( ) const
pure virtual

◆ ConstNodeBegin() [1/3]

virtual ConstNodeIterator VAPoR::Grid::ConstNodeBegin ( ) const
inlinevirtual

Return constant grid node coordinate iterator

If minu and maxu are specified the iterator is constrained to operation within the axis-aligned box defined by minu and maxu.

Parameters
[in]minuMinimum box coordinate.
[in]maxuMaximum box coordinate.

Definition at line 1021 of file Grid.h.

◆ ConstNodeBegin() [2/3]

virtual ConstNodeIterator VAPoR::Grid::ConstNodeBegin ( const CoordType minu,
const CoordType maxu 
) const
inlinevirtual

Definition at line 1023 of file Grid.h.

◆ ConstNodeBegin() [3/3]

virtual ConstNodeIterator VAPoR::Grid::ConstNodeBegin ( const std::vector< double > &  minu,
const std::vector< double > &  maxu 
) const
inlinevirtual

Definition at line 1028 of file Grid.h.

◆ ConstNodeEnd()

virtual ConstNodeIterator VAPoR::Grid::ConstNodeEnd ( ) const
inlinevirtual

Definition at line 1036 of file Grid.h.

◆ CopyFromArr3() [1/2]

template<typename T >
static void VAPoR::Grid::CopyFromArr3 ( const std::array< T, 3 > &  src,
std::vector< T > &  dst 
)
inlinestatic

Definition at line 1225 of file Grid.h.

◆ CopyFromArr3() [2/2]

template<typename T >
static void VAPoR::Grid::CopyFromArr3 ( const std::array< T, 3 > &  src,
T *  dst 
)
inlinestatic

Definition at line 1230 of file Grid.h.

◆ CopyToArr3() [1/2]

template<typename T >
static void VAPoR::Grid::CopyToArr3 ( const std::vector< T > &  src,
std::array< T, 3 > &  dst 
)
inlinestatic

Definition at line 1216 of file Grid.h.

◆ CopyToArr3() [2/2]

template<typename T >
static void VAPoR::Grid::CopyToArr3 ( const T *  src,
size_t  n,
std::array< T, 3 > &  dst 
)
inlinestatic

Definition at line 1221 of file Grid.h.

◆ Dims()

static DimsType VAPoR::Grid::Dims ( const DimsType min,
const DimsType max 
)
static

Compute the dimensions of a rectangular region bounded by min and max coordinates

◆ end()

Iterator VAPoR::Grid::end ( )
inline

Definition at line 1201 of file Grid.h.

◆ GetBlks()

const std::vector< float * > & VAPoR::Grid::GetBlks ( ) const
inline

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

Definition at line 168 of file Grid.h.

◆ GetBlockSize()

const std::vector< size_t > & VAPoR::Grid::GetBlockSize ( ) const
inline

Return the internal blocking factor

This method returns the internal blocking factor passed to the constructor.

Definition at line 163 of file Grid.h.

◆ GetBoundingBox() [1/2]

virtual void VAPoR::Grid::GetBoundingBox ( const DimsType min,
const DimsType max,
CoordType minu,
CoordType maxu 
) const
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.

Note
The results are undefined if any index of 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()

Parameters
[in]minAn array containing the minimum grid indices (offsets).
[in]maxAn array containing the maximum grid indices (offsets).
[out]minuA two-to-three element array containing the minimum user coordinates.
[out]maxuA two-to-three element array containing the maximum user coordinates.
See also
GetDimensions(), Grid()

Implemented in VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.

◆ GetBoundingBox() [2/2]

virtual void VAPoR::Grid::GetBoundingBox ( const std::vector< size_t > &  min,
const std::vector< size_t > &  max,
std::vector< double > &  minu,
std::vector< double > &  maxu 
) const
inlinevirtual

Definition at line 348 of file Grid.h.

References VAssert.

◆ GetCellDimensions()

virtual const DimsType & VAPoR::Grid::GetCellDimensions ( ) const
pure virtual

◆ GetCellNeighbors()

virtual bool VAPoR::Grid::GetCellNeighbors ( const DimsType cindices,
std::vector< DimsType > &  cells 
) const
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

Parameters
[in]cindicesAn ordered vector of multi-dimensional cell indices.
[out]cellsA vector of index vectors. Each index vector has size given by GetDimensions.size()

Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.

◆ GetCellNodes() [1/2]

virtual bool VAPoR::Grid::GetCellNodes ( const DimsType cindices,
std::vector< DimsType > &  nodes 
) const
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.

Parameters
[in]cindicesAn array of size Grid::GetDimensions.size() specifying the index of the cell to query.
[out]nodesAn array of size Grid::GetMaxVertexPerCell() x Grid::GetDimensions.size() that will contain the concatenated list of node indices on return.
[out]nThe number of node indices returned in nodes, an integer in the range (0..Grid::GetMaxVertexPerCell()).

Implemented in VAPoR::StructuredGrid, VAPoR::UnstructuredGrid, and VAPoR::StructuredGrid.

◆ GetCellNodes() [2/2]

virtual bool VAPoR::Grid::GetCellNodes ( const size_t  cindices[],
std::vector< DimsType > &  nodes 
) const
inlinevirtual
Deprecated:

Reimplemented in VAPoR::StructuredGrid.

Definition at line 609 of file Grid.h.

◆ GetCellOffset()

virtual long VAPoR::Grid::GetCellOffset ( ) const
inlinevirtual

Get the linear offset to the cell IDs

Return the smallest Cell ID. The default is zero

Definition at line 742 of file Grid.h.

◆ GetCoordDimensions()

virtual DimsType VAPoR::Grid::GetCoordDimensions ( size_t  dim) const
pure virtual

Return the dimensions of the specified coordinate variable

Parameters
[in]dimAn 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.

◆ GetDimensionInBlks()

const std::vector< size_t > VAPoR::Grid::GetDimensionInBlks ( ) const
inline

Return the ijk dimensions of grid in blocks

Returns the number of blocks defined along each axis of the grid

Parameters
[out]dimsA two or three element array containing the grid dimension in blocks
See also
GetBlockSize();

Definition at line 156 of file Grid.h.

◆ GetDimensions()

const DimsType & VAPoR::Grid::GetDimensions ( ) const
inline

Return the dimensions of grid connectivity array

Parameters
[out]dimsThe value of dims parameter provided to the constructor. If the parameter has less than 3 values, then number 1 will be filled.
See also
GetGeometryDim(), GetTopologyDim()

Definition at line 100 of file Grid.h.

◆ GetEnclosingRegion()

virtual bool VAPoR::Grid::GetEnclosingRegion ( const CoordType minu,
const CoordType maxu,
DimsType min,
DimsType max 
) const
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

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

Implemented in VAPoR::StructuredGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.

◆ GetGeometryDim()

virtual size_t VAPoR::Grid::GetGeometryDim ( ) const
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.

See also
GetTopologyDim()

Implemented in VAPoR::CurvilinearGrid, VAPoR::LayeredGrid, VAPoR::RegularGrid, VAPoR::StretchedGrid, VAPoR::UnstructuredGrid2D, VAPoR::UnstructuredGrid3D, VAPoR::UnstructuredGridCoordless, and VAPoR::UnstructuredGridLayered.

◆ GetIndicesCell() [1/3]

virtual bool VAPoR::Grid::GetIndicesCell ( const CoordType coords,
DimsType indices 
) const
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.

Return values
statustrue 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.

◆ GetIndicesCell() [2/3]

virtual bool VAPoR::Grid::GetIndicesCell ( const double  coords[3],
size_t  indices[3] 
) const
inlinevirtual
Deprecated:

Reimplemented in VAPoR::LayeredGrid, and VAPoR::UnstructuredGrid2D.

Definition at line 511 of file Grid.h.

◆ GetIndicesCell() [3/3]

virtual bool VAPoR::Grid::GetIndicesCell ( const std::vector< double > &  coords,
std::vector< size_t > &  indices 
) const
inlinevirtual
Deprecated:

Reimplemented in VAPoR::LayeredGrid, and VAPoR::UnstructuredGrid2D.

Definition at line 523 of file Grid.h.

◆ GetInterpolationOrder()

virtual int VAPoR::Grid::GetInterpolationOrder ( ) const
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

See also
SetInterpolationOrder()

Reimplemented in VAPoR::LayeredGrid.

Definition at line 431 of file Grid.h.

◆ GetMaxVertexPerCell()

virtual size_t VAPoR::Grid::GetMaxVertexPerCell ( ) const
pure virtual

Return the maximum number of vertices per cell

Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.

◆ GetMaxVertexPerFace()

virtual size_t VAPoR::Grid::GetMaxVertexPerFace ( ) const
pure virtual

Return the maximum number of vertices per cell face

Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.

◆ GetMinAbs()

virtual DimsType VAPoR::Grid::GetMinAbs ( ) const
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.

Note
The value of returned is not used within the Grid class and any value can be stored here using SetMinAbs().

Definition at line 755 of file Grid.h.

◆ GetMissingValue()

virtual float VAPoR::Grid::GetMissingValue ( ) const
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.

◆ GetNodeCells()

virtual bool VAPoR::Grid::GetNodeCells ( const DimsType indices,
std::vector< DimsType > &  cells 
) const
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

Parameters
[out]nodesA vector of index vectors . Each index vector has size given by GetDimensions.size()

Implemented in VAPoR::StructuredGrid, and VAPoR::UnstructuredGrid.

◆ GetNodeDimensions()

virtual const DimsType & VAPoR::Grid::GetNodeDimensions ( ) const
pure virtual

◆ GetNodeOffset()

virtual long VAPoR::Grid::GetNodeOffset ( ) const
inlinevirtual

Get the linear offset to the node IDs

Return the smallest node ID. The default is zero

Definition at line 735 of file Grid.h.

◆ GetNumCellDimensions()

virtual const size_t VAPoR::Grid::GetNumCellDimensions ( ) const
pure virtual

◆ GetNumDimensions() [1/2]

size_t VAPoR::Grid::GetNumDimensions ( ) const
inline

Return the useful number of dimensions of grid connectivity array

Parameters
[out]dimsThe 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().

◆ GetNumDimensions() [2/2]

static size_t VAPoR::Grid::GetNumDimensions ( DimsType  dims)
static

Return the number of non-unit dimensions

Return the number of non-unit dimensions in dims

Parameters
[in]dims

◆ GetNumNodeDimensions()

virtual const size_t VAPoR::Grid::GetNumNodeDimensions ( ) const
pure virtual

◆ GetPeriodic()

virtual const std::vector< bool > & VAPoR::Grid::GetPeriodic ( ) const
inlinevirtual

Check for periodic boundaries

Definition at line 718 of file Grid.h.

◆ GetRange() [1/3]

virtual void VAPoR::Grid::GetRange ( const DimsType min,
const DimsType max,
float  range[2] 
) const
virtual

◆ GetRange() [2/3]

virtual void VAPoR::Grid::GetRange ( float  range[2]) const
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.

Parameters
[out]range[2]A two-element array containing the mininum and maximum values, in that order

◆ GetRange() [3/3]

virtual void VAPoR::Grid::GetRange ( std::vector< size_t >  min,
std::vector< size_t >  max,
float  range[2] 
) const
inlinevirtual
Deprecated:

Definition at line 549 of file Grid.h.

◆ GetTopologyDim()

virtual size_t VAPoR::Grid::GetTopologyDim ( ) const
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.

See also
GetGeometryDim()

Reimplemented in VAPoR::ConstantGrid.

Definition at line 729 of file Grid.h.

◆ GetType()

virtual std::string VAPoR::Grid::GetType ( ) const
pure virtual

◆ GetUserCoordinates() [1/6]

virtual void VAPoR::Grid::GetUserCoordinates ( const DimsType indices,
CoordType coords 
) const
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.

Parameters
[in]indexLogical 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]coordUser 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.

◆ GetUserCoordinates() [2/6]

virtual void VAPoR::Grid::GetUserCoordinates ( const size_t  indices[],
double  coords[] 
) const
inlinevirtual

Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.

Definition at line 469 of file Grid.h.

◆ GetUserCoordinates() [3/6]

virtual void VAPoR::Grid::GetUserCoordinates ( const std::vector< size_t > &  indices,
std::vector< double > &  coords 
) const
inlinevirtual

Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.

Definition at line 478 of file Grid.h.

◆ GetUserCoordinates() [4/6]

virtual void VAPoR::Grid::GetUserCoordinates ( size_t  i,
double &  x,
double &  y,
double &  z 
) const
virtual

Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.

◆ GetUserCoordinates() [5/6]

virtual void VAPoR::Grid::GetUserCoordinates ( size_t  i,
size_t  j,
double &  x,
double &  y,
double &  z 
) const
virtual

Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.

◆ GetUserCoordinates() [6/6]

virtual void VAPoR::Grid::GetUserCoordinates ( size_t  i,
size_t  j,
size_t  k,
double &  x,
double &  y,
double &  z 
) const
virtual

Reimplemented in VAPoR::LayeredGrid, and VAPoR::RegularGrid.

◆ GetUserExtents() [1/3]

virtual void VAPoR::Grid::GetUserExtents ( CoordType minu,
CoordType maxu 
) const
virtual

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.

Parameters
[out]minuA two or three element array containing the minimum user coordinates.
[out]maxuA two or three element array containing the maximum user coordinates.
See also
GetDimensions(), Grid()

◆ GetUserExtents() [2/3]

virtual void VAPoR::Grid::GetUserExtents ( double  minu[3],
double  maxu[3] 
) const
inlinevirtual
Deprecated:

Definition at line 293 of file Grid.h.

◆ GetUserExtents() [3/3]

virtual void VAPoR::Grid::GetUserExtents ( std::vector< double > &  minu,
std::vector< double > &  maxu 
) const
inlinevirtual
Deprecated:

Definition at line 304 of file Grid.h.

◆ GetUserExtentsHelper()

virtual void VAPoR::Grid::GetUserExtentsHelper ( CoordType minu,
CoordType maxu 
) const
protectedpure virtual

◆ GetValue() [1/5]

virtual float VAPoR::Grid::GetValue ( const CoordType coords) 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]coordsA vector of size matching the topology dimension of the mesh whose contents specify the coordinates of a point in space.
See also
GetInterpolationOrder(), HasPeriodic(), GetMissingValue()
GetUserExtents()

Reimplemented in VAPoR::ConstantGrid, and VAPoR::LayeredGrid.

◆ GetValue() [2/5]

virtual float VAPoR::Grid::GetValue ( const double  coords[]) const
inlinevirtual
Deprecated:

Definition at line 257 of file Grid.h.

◆ GetValue() [3/5]

virtual float VAPoR::Grid::GetValue ( const std::vector< double > &  coords) const
inlinevirtual
Deprecated:

Definition at line 248 of file Grid.h.

◆ GetValue() [4/5]

virtual float VAPoR::Grid::GetValue ( double  x,
double  y 
) const
inlinevirtual

Definition at line 264 of file Grid.h.

◆ GetValue() [5/5]

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

Reimplemented in VAPoR::SphericalGrid.

Definition at line 269 of file Grid.h.

◆ GetValueAtIndex() [1/2]

virtual float VAPoR::Grid::GetValueAtIndex ( const DimsType indices) const
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

Parameters
[in]indicesof grid point along fastest varying dimension.

◆ GetValueAtIndex() [2/2]

virtual float VAPoR::Grid::GetValueAtIndex ( const std::vector< size_t > &  indices) const
inlinevirtual
Deprecated:

Definition at line 186 of file Grid.h.

◆ GetValueLinear()

virtual float VAPoR::Grid::GetValueLinear ( const CoordType coords) const
protectedpure virtual

◆ GetValueNearestNeighbor()

virtual float VAPoR::Grid::GetValueNearestNeighbor ( const CoordType coords) const
protectedpure virtual

◆ GetValuePtrAtIndex()

virtual float * VAPoR::Grid::GetValuePtrAtIndex ( const std::vector< float * > &  blks,
const DimsType indices 
) const
protectedvirtual

◆ HasInvertedCoordinateSystemHandiness()

virtual bool VAPoR::Grid::HasInvertedCoordinateSystemHandiness ( ) const
inlinevirtual

Return true if mesh primitives have counter clockwise winding order.

Reimplemented in VAPoR::StructuredGrid.

Definition at line 422 of file Grid.h.

◆ HasMissingData()

bool VAPoR::Grid::HasMissingData ( ) const
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.

Definition at line 417 of file Grid.h.

◆ InsideGrid() [1/3]

virtual bool VAPoR::Grid::InsideGrid ( const CoordType coords) const
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

Parameters
[in]coordsUser coordinates of point in space
Return values
boolTrue 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.

◆ InsideGrid() [2/3]

virtual bool VAPoR::Grid::InsideGrid ( const double  coords[3]) const
inlinevirtual
Deprecated:

Definition at line 571 of file Grid.h.

◆ InsideGrid() [3/3]

virtual bool VAPoR::Grid::InsideGrid ( const std::vector< double > &  coords) const
inlinevirtual
Deprecated:

Definition at line 580 of file Grid.h.

◆ PointInsideBoundingRectangle()

static bool VAPoR::Grid::PointInsideBoundingRectangle ( const double  pt[],
const double  verts[],
int  n 
)
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

Parameters
[in]ptA two-element array of point coordinates
[in]vertsAn array of dimensions n * 2 containing a list of points.
[in]nThe number of 2D points in verts
Return values
statusReturns true if \pt is inside or on the bounding rectangle of verts. False otherwise.

Definition at line 781 of file Grid.h.

References VAssert.

◆ SetCellOffset()

virtual void VAPoR::Grid::SetCellOffset ( long  offset)
inlinevirtual

Reimplemented in VAPoR::UnstructuredGrid.

Definition at line 743 of file Grid.h.

◆ SetHasMissingValues()

void VAPoR::Grid::SetHasMissingValues ( bool  flag)
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.

Definition at line 407 of file Grid.h.

◆ SetInterpolationOrder()

virtual void VAPoR::Grid::SetInterpolationOrder ( int  order)
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

Parameters
[in]orderinterpolation order
See also
GetInterpolationOrder()

Reimplemented in VAPoR::LayeredGrid.

◆ SetMinAbs()

virtual void VAPoR::Grid::SetMinAbs ( const DimsType minAbs)
inlinevirtual

Set the absolute minimum grid coordinate

Parameters
[in]minAbsMust 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

Definition at line 764 of file Grid.h.

◆ SetMissingValue()

void VAPoR::Grid::SetMissingValue ( float  missing_value)
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.

See also
HasMissingData(), GetMissingValue()

Definition at line 398 of file Grid.h.

◆ SetNodeOffset()

virtual void VAPoR::Grid::SetNodeOffset ( long  offset)
inlinevirtual

Reimplemented in VAPoR::UnstructuredGrid.

Definition at line 736 of file Grid.h.

◆ SetPeriodic()

virtual void VAPoR::Grid::SetPeriodic ( const std::vector< bool > &  periodic)
inlinevirtual

Set periodic boundaries

This method changes the periodicity of boundaries set by the class constructor

Parameters
[in]periodicA boolean vector of size given by GetGeometryDim() indicating which coordinate axes are periodic.

Reimplemented in VAPoR::LayeredGrid.

Definition at line 707 of file Grid.h.

◆ SetValue() [1/2]

virtual void VAPoR::Grid::SetValue ( const DimsType indices,
float  value 
)
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.

◆ SetValue() [2/2]

virtual void VAPoR::Grid::SetValue ( const size_t  indices[3],
float  value 
)
inlinevirtual
Deprecated:

Definition at line 202 of file Grid.h.

◆ SetValueIJK() [1/3]

void VAPoR::Grid::SetValueIJK ( size_t  i,
float  v 
)
inline

Definition at line 222 of file Grid.h.

References SetValueIJK().

Referenced by SetValueIJK().

◆ SetValueIJK() [2/3]

void VAPoR::Grid::SetValueIJK ( size_t  i,
size_t  j,
float  v 
)
inline

Definition at line 221 of file Grid.h.

References SetValueIJK().

Referenced by SetValueIJK().

◆ SetValueIJK() [3/3]

void VAPoR::Grid::SetValueIJK ( size_t  i,
size_t  j,
size_t  k,
float  v 
)

◆ TrilinearInterpolate()

float VAPoR::Grid::TrilinearInterpolate ( size_t  i,
size_t  j,
size_t  k,
const double  xwgt,
const double  ywgt,
const double  zwgt 
) const
protected

Friends And Related Function Documentation

◆ operator<<

VDF_API friend std::ostream & operator<< ( std::ostream &  o,
const Grid g 
)
friend

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