VAPOR3 3.9.4
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | Friends | List of all members
VAPoR::UnstructuredGrid Class Reference

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

#include <UnstructuredGrid.h>

Inheritance diagram for VAPoR::UnstructuredGrid:
VAPoR::Grid VAPoR::UnstructuredGrid2D VAPoR::UnstructuredGrid3D VAPoR::UnstructuredGridCoordless VAPoR::UnstructuredGridLayered

Public Types

enum  Location { NODE , CELL , EDGE }
 
- Public Types inherited from VAPoR::Grid
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

 UnstructuredGrid (const DimsType &vertexDims, const DimsType &faceDims, const DimsType &edgeDims, const DimsType &bs, const std::vector< float * > &blks, size_t topology_dimension, const int *vertexOnFace, const int *faceOnVertex, const int *faceOnFace, Location location, size_t maxVertexPerFace, size_t maxFacePerVertex, long nodeOffset, long cellOffset)
 
 UnstructuredGrid (const std::vector< size_t > &vertexDims, const std::vector< size_t > &faceDims, const std::vector< size_t > &edgeDims, const std::vector< size_t > &bs, const std::vector< float * > &blks, size_t topology_dimension, const int *vertexOnFace, const int *faceOnVertex, const int *faceOnFace, Location location, size_t maxVertexPerFace, size_t maxFacePerVertex, long nodeOffset, long cellOffset)
 
 UnstructuredGrid ()=default
 
virtual ~UnstructuredGrid ()=default
 
std::string GetType () const override
 
bool GetCellNodes (const DimsType &cindices, std::vector< DimsType > &nodes) const override
 
virtual bool GetCellNeighbors (const DimsType &cindices, std::vector< DimsType > &cells) const override
 
virtual bool GetNodeCells (const DimsType &indices, std::vector< DimsType > &cells) const override
 
size_t GetMaxVertexPerFace () const override
 
size_t GetMaxVertexPerCell () const override
 
const VAPoR::DimsTypeGetNodeDimensions () const override
 
const size_t GetNumNodeDimensions () const override
 
const DimsTypeGetCellDimensions () const override
 
virtual const size_t GetNumCellDimensions () const override
 
const DimsTypeGetEdgeDimensions () const
 
size_t GetMissingID () const
 
void SetMissingID (size_t v)
 
size_t GetBoundaryID () const
 
void SetBoundaryID (size_t v)
 
virtual void ClampCoord (const CoordType &coords, CoordType &cCoords) const override
 
virtual void ClampCoord (const double coords[3], double cCoords[3]) const override
 
virtual void SetNodeOffset (long offset) override
 
virtual void SetCellOffset (long offset) override
 
- Public Member Functions inherited from VAPoR::Grid
 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 std::string GetClassType ()
 
- Static Public Member Functions inherited from VAPoR::Grid
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 Attributes

const int * _vertexOnFace
 
const int * _faceOnVertex
 
const int * _faceOnFace
 
size_t _maxVertexPerFace
 
size_t _maxFacePerVertex
 
Location _location
 

Friends

VDF_API friend std::ostream & operator<< (std::ostream &o, const UnstructuredGrid &sg)
 

Additional Inherited Members

- Protected Member Functions inherited from VAPoR::Grid
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
 

Detailed Description

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

Author
John Clyne

This abstract base class defines a 2D or 3D unstructured grid.

The unstructured grid samples a scalar function at each grid point. The scalar function samples are stored as an array decomposed into equal-sized blocks.

Because grid samples are repesented internally as arrays, when accessing multiple grid points better performance is achieved by using unit stride. The I axis varies fastest (has unit stride), followed by J, then K. Best performance is achieved when using the class iterator: UnstructuredGrid::Iterator.

For methods that allow the specification of grid indecies or coordinates as a single parameter tuple (e.g. vector <double> coordinate) the first element of the tuple corresponds to the I axis, etc.

Note
Throughout this class grid vertex offsets are specified as i, j, k, where i, j, k are integers. User coordinates are real values denoted x, y, z, and are given by functions X(i,j,k), Y(i,j,k), Z(i,j,k).

Definition at line 43 of file UnstructuredGrid.h.

Member Enumeration Documentation

◆ Location

Enumerator
NODE 
CELL 
EDGE 

Definition at line 45 of file UnstructuredGrid.h.

Constructor & Destructor Documentation

◆ UnstructuredGrid() [1/3]

VAPoR::UnstructuredGrid::UnstructuredGrid ( const DimsType vertexDims,
const DimsType faceDims,
const DimsType edgeDims,
const DimsType bs,
const std::vector< float * > &  blks,
size_t  topology_dimension,
const int *  vertexOnFace,
const int *  faceOnVertex,
const int *  faceOnFace,
Location  location,
size_t  maxVertexPerFace,
size_t  maxFacePerVertex,
long  nodeOffset,
long  cellOffset 
)

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

Parameters
[in]vertexDimsDimensions of grid nodes (vertices). The product of the elements of vertexDims gives the total number of nodes in the mesh. If the rank of vertexDims is greater than one the mesh is assumed to be structured along the slower varying dimensions (e.g. a layered mesh)
[in]faceDimsDimensions of grid cells. The product of the elements of cells gives the total number of cells in the mesh. If the rank of cells is greater than one the mesh is assumed to be structured along the slower varying dimensions (e.g. a layered mesh). Moreover, the slowest varying dimensions must be one less than the corresponding dimension of vertexDims.
[in]edgeDimsDimensions of grid edges. The product of the elements of edges gives the total number of edges in the mesh. If the rank of edges is greater than one the mesh is assumed to be structured along the slower varying dimensions (e.g. a layered mesh). Moreover, the slowest varying dimensions must be one less than the corresponding dimension of vertexDims.
[in]blksGrid data values. The location of the data values (node, cell, edge) is determined by location. If the dimensions (vertexDims, etc) are multi-dimensional the size of blks must match the size of the slowest varying dimension. Each element of blks must point to an area of memory of size n elements, where n is the first element of vertexDims, faceDims, or edgeDims as indicated by location.
[in]vertexOnFaceAn array with dimensions:

faceDims[0] * maxVertexPerFace that provides for each cell the 1D node IDs of each corner node. If the number of corner nodes is less than maxVertexPerFace the missing node indices will be equal to GetMissingID(); The ordering of the nodes is counter-clockwise.

Parameters
[in]faceOnVertexAn array with dimensions:

vertexDims[0] * maxFacePerVertex that provides for each vertex the 1D node IDs of each face sharing that vertex. If the number of faces is less than maxFacePerVertex the missing face indices will be equal to GetMissingID(); The ordering of the faces is counter-clockwise.

Parameters
[in]faceOnFaceAn array with dimensions:

faceDims[0] * maxVertexPerFace that provides for each cell the 1D cell IDs of border cell If the number of corner nodes is less than maxVertexPerFace the missing node indices will be equal to GetMissingID(). If an edge is on the mesh boundary the index will be set to the value GetBoundaryIndex(). The ordering of the neighboring faces is counter-clockwise.

Parameters
[in]locationThe location of grid data: at the nodes, edge-centered, or cell-centered
[in]maxVertexPerFaceThe maxium number of nodes that a face may have.
[in]nodeOffsetThe offset from zero for the first element in vertexOnFace
[in]cellOffsetThe offset from zero for the first element in faceOnVertex or faceOnFace

◆ UnstructuredGrid() [2/3]

VAPoR::UnstructuredGrid::UnstructuredGrid ( const std::vector< size_t > &  vertexDims,
const std::vector< size_t > &  faceDims,
const std::vector< size_t > &  edgeDims,
const std::vector< size_t > &  bs,
const std::vector< float * > &  blks,
size_t  topology_dimension,
const int *  vertexOnFace,
const int *  faceOnVertex,
const int *  faceOnFace,
Location  location,
size_t  maxVertexPerFace,
size_t  maxFacePerVertex,
long  nodeOffset,
long  cellOffset 
)

◆ UnstructuredGrid() [3/3]

VAPoR::UnstructuredGrid::UnstructuredGrid ( )
default

◆ ~UnstructuredGrid()

virtual VAPoR::UnstructuredGrid::~UnstructuredGrid ( )
virtualdefault

Member Function Documentation

◆ ClampCoord() [1/2]

virtual void VAPoR::UnstructuredGrid::ClampCoord ( const CoordType coords,
CoordType cCoords 
) const
inlineoverridevirtual

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()

Implements VAPoR::Grid.

Definition at line 179 of file UnstructuredGrid.h.

◆ ClampCoord() [2/2]

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

Reimplemented from VAPoR::Grid.

Definition at line 183 of file UnstructuredGrid.h.

◆ GetBoundaryID()

size_t VAPoR::UnstructuredGrid::GetBoundaryID ( ) const
inline

Get boundary element ID

Return the value used to indicate termination of a list of element IDs

Definition at line 176 of file UnstructuredGrid.h.

◆ GetCellDimensions()

const DimsType & VAPoR::UnstructuredGrid::GetCellDimensions ( ) const
inlineoverridevirtual

Return the grid cell dimmensions

Implements VAPoR::Grid.

Definition at line 157 of file UnstructuredGrid.h.

◆ GetCellNeighbors()

virtual bool VAPoR::UnstructuredGrid::GetCellNeighbors ( const DimsType cindices,
std::vector< DimsType > &  cells 
) const
overridevirtual

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()

Implements VAPoR::Grid.

◆ GetCellNodes()

bool VAPoR::UnstructuredGrid::GetCellNodes ( const DimsType cindices,
std::vector< DimsType > &  nodes 
) const
overridevirtual

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()).

Implements VAPoR::Grid.

◆ GetClassType()

static std::string VAPoR::UnstructuredGrid::GetClassType ( )
inlinestatic

Definition at line 133 of file UnstructuredGrid.h.

◆ GetEdgeDimensions()

const DimsType & VAPoR::UnstructuredGrid::GetEdgeDimensions ( ) const
inline

Return the grid edge dimmensions

Definition at line 163 of file UnstructuredGrid.h.

◆ GetMaxVertexPerCell()

size_t VAPoR::UnstructuredGrid::GetMaxVertexPerCell ( ) const
inlineoverridevirtual

Return the maximum number of vertices per cell

Implements VAPoR::Grid.

Definition at line 148 of file UnstructuredGrid.h.

◆ GetMaxVertexPerFace()

size_t VAPoR::UnstructuredGrid::GetMaxVertexPerFace ( ) const
inlineoverridevirtual

Return the maximum number of vertices per cell face

Implements VAPoR::Grid.

Definition at line 146 of file UnstructuredGrid.h.

◆ GetMissingID()

size_t VAPoR::UnstructuredGrid::GetMissingID ( ) const
inline

Get missing element ID

Return the value used to indicate termination of a list of element IDs

Definition at line 169 of file UnstructuredGrid.h.

◆ GetNodeCells()

virtual bool VAPoR::UnstructuredGrid::GetNodeCells ( const DimsType indices,
std::vector< DimsType > &  cells 
) const
overridevirtual

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()

Implements VAPoR::Grid.

◆ GetNodeDimensions()

const VAPoR::DimsType & VAPoR::UnstructuredGrid::GetNodeDimensions ( ) const
overridevirtual

Return the grid node dimmensions

Implements VAPoR::Grid.

◆ GetNumCellDimensions()

virtual const size_t VAPoR::UnstructuredGrid::GetNumCellDimensions ( ) const
inlineoverridevirtual

Implements VAPoR::Grid.

Definition at line 159 of file UnstructuredGrid.h.

◆ GetNumNodeDimensions()

const size_t VAPoR::UnstructuredGrid::GetNumNodeDimensions ( ) const
overridevirtual

Implements VAPoR::Grid.

◆ GetType()

std::string VAPoR::UnstructuredGrid::GetType ( ) const
inlineoverridevirtual

◆ SetBoundaryID()

void VAPoR::UnstructuredGrid::SetBoundaryID ( size_t  v)
inline

Definition at line 177 of file UnstructuredGrid.h.

◆ SetCellOffset()

virtual void VAPoR::UnstructuredGrid::SetCellOffset ( long  offset)
inlineoverridevirtual

Reimplemented from VAPoR::Grid.

Definition at line 188 of file UnstructuredGrid.h.

◆ SetMissingID()

void VAPoR::UnstructuredGrid::SetMissingID ( size_t  v)
inline

Definition at line 170 of file UnstructuredGrid.h.

◆ SetNodeOffset()

virtual void VAPoR::UnstructuredGrid::SetNodeOffset ( long  offset)
inlineoverridevirtual

Reimplemented from VAPoR::Grid.

Definition at line 187 of file UnstructuredGrid.h.

Friends And Related Function Documentation

◆ operator<<

VDF_API friend std::ostream & operator<< ( std::ostream &  o,
const UnstructuredGrid sg 
)
friend

Member Data Documentation

◆ _faceOnFace

const int* VAPoR::UnstructuredGrid::_faceOnFace
protected

Definition at line 201 of file UnstructuredGrid.h.

◆ _faceOnVertex

const int* VAPoR::UnstructuredGrid::_faceOnVertex
protected

Definition at line 200 of file UnstructuredGrid.h.

◆ _location

Location VAPoR::UnstructuredGrid::_location
protected

Definition at line 204 of file UnstructuredGrid.h.

◆ _maxFacePerVertex

size_t VAPoR::UnstructuredGrid::_maxFacePerVertex
protected

Definition at line 203 of file UnstructuredGrid.h.

◆ _maxVertexPerFace

size_t VAPoR::UnstructuredGrid::_maxVertexPerFace
protected

Definition at line 202 of file UnstructuredGrid.h.

◆ _vertexOnFace

const int* VAPoR::UnstructuredGrid::_vertexOnFace
protected

Definition at line 199 of file UnstructuredGrid.h.


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