VAPOR3 3.9.4
Public Member Functions | Static Public Member Functions | List of all members
VAPoR::SphericalGrid Class Reference

This class implements a 2D or 3D spherical grid. More...

#include <SphericalGrid.h>

Inheritance diagram for VAPoR::SphericalGrid:
VAPoR::RegularGrid VAPoR::StructuredGrid VAPoR::Grid

Public Member Functions

 SphericalGrid (const size_t bs[3], const size_t min[3], const size_t max[3], const double extents[6], const size_t permutation[3], const bool periodic[3], float **blks)
 
 SphericalGrid (const size_t bs[3], const size_t min[3], const size_t max[3], const double extents[6], const size_t permutation[3], const bool periodic[3], float **blks, float missing_value)
 
float GetValue (double x, double y, double z) const
 
virtual void GetUserExtents (double extents[6]) const
 
virtual void GetBoundingBox (const size_t min[3], const size_t max[3], double extents[6]) const
 
int GetUserCoordinates (size_t i, size_t j, size_t k, double *x, double *y, double *z) const
 
void GetIJKIndex (double x, double y, double z, size_t *i, size_t *j, size_t *k) const
 
void GetIJKIndexFloor (double x, double y, double z, size_t *i, size_t *j, size_t *k) const
 
bool InsideGrid (double x, double y, double z) const
 
- Public Member Functions inherited from VAPoR::RegularGrid
 RegularGrid (const DimsType &dims, const DimsType &bs, const std::vector< float * > &blks, const CoordType &minu, const CoordType &maxu)
 
 RegularGrid (const std::vector< size_t > &dims, const std::vector< size_t > &bs, const std::vector< float * > &blks, const std::vector< double > &minu, const std::vector< double > &maxu)
 
 RegularGrid ()=default
 
virtual ~RegularGrid ()=default
 
virtual size_t GetGeometryDim () const override
 
virtual DimsType GetCoordDimensions (size_t dim) const override
 
std::string GetType () const override
 
virtual void GetBoundingBox (const DimsType &min, const DimsType &max, CoordType &minu, CoordType &maxu) const override
 
virtual void GetUserCoordinates (const DimsType &indices, CoordType &coords) const override
 
virtual bool GetIndicesCell (const CoordType &coords, DimsType &indices) const override
 
virtual bool InsideGrid (const CoordType &coords) const override
 
virtual ConstCoordItr ConstCoordBegin () const override
 Return constant grid coordinate iterator.
 
virtual ConstCoordItr ConstCoordEnd () const override
 
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
 
- Public Member Functions inherited from VAPoR::StructuredGrid
 StructuredGrid (const std::vector< size_t > &dims, const std::vector< size_t > &bs, const std::vector< float * > &blks)
 
 StructuredGrid (const DimsType &dims, const DimsType &bs, const std::vector< float * > &blks)
 
 StructuredGrid ()=default
 
virtual ~StructuredGrid ()=default
 
std::string GetType () const override
 
const DimsTypeGetNodeDimensions () const override
 
const size_t GetNumNodeDimensions () const override
 
const DimsTypeGetCellDimensions () const override
 
const size_t GetNumCellDimensions () const override
 
virtual 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 &cindices, std::vector< DimsType > &cells) const override
 
virtual bool GetEnclosingRegion (const CoordType &minu, const CoordType &maxu, DimsType &min, DimsType &max) const override
 
size_t GetMaxVertexPerFace () const override
 
size_t GetMaxVertexPerCell () const override
 
virtual void ClampCoord (const CoordType &coords, CoordType &cCoords) const override
 
virtual void ClampCoord (const double coords[3], double cCoords[3]) const override
 
virtual bool HasInvertedCoordinateSystemHandiness () const override
 
virtual bool GetCellNodes (const DimsType &cindices, std::vector< DimsType > &nodes) const=0
 
virtual bool GetCellNodes (const size_t cindices[], std::vector< DimsType > &nodes) const
 
- 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 void CartToSph (double x, double y, double z, double *phi, double *theta, double *r)
 
static void SphToCart (double phi, double theta, double r, double *x, double *y, double *z)
 
- Static Public Member Functions inherited from VAPoR::RegularGrid
static std::string GetClassType ()
 
- Static Public Member Functions inherited from VAPoR::StructuredGrid
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)
 

Additional Inherited Members

- 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
 
- Protected Member Functions inherited from VAPoR::RegularGrid
virtual float GetValueNearestNeighbor (const CoordType &coords) const override
 
virtual float GetValueLinear (const CoordType &coords) const override
 
virtual void GetUserExtentsHelper (CoordType &minu, CoordType &maxu) const override
 
- 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

This class implements a 2D or 3D spherical grid.

This class implements a 2D or 3D spherical grid: a generalization of a regular grid where the spacing of grid points along a single dimension may vary at each grid point. The spacing along the remaining one (2D case) or two (3D case) dimensions is invariant between grid points. For example, if K is the spherical dimension than the z coordinate is given by some function f(i,j,k):

z = f(i,j,k)

while the remaining x and y coordinates are givey by (i*dx, j*dy) for some real dx and dy .

Definition at line 29 of file SphericalGrid.h.

Constructor & Destructor Documentation

◆ SphericalGrid() [1/2]

VAPoR::SphericalGrid::SphericalGrid ( const size_t  bs[3],
const size_t  min[3],
const size_t  max[3],
const double  extents[6],
const size_t  permutation[3],
const bool  periodic[3],
float **  blks 
)

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

Parameters
[in]bsA three-element vector specifying the dimensions of each block storing the sampled scalar function.
[in]minA three-element vector specifying the ijk index of the first point in the grid. The first grid point need not coincide with block boundaries. I.e. the indecies need not be (0,0,0)
[in]maxA three-element vector specifying the ijk index of the last point in the grid
[in]extentsA six-element vector specifying the user coordinates of the first (first three elements) and last (last three elements) of the grid points indicated by min and max, respectively. The units are degrees, not radians.
[in]permutationA three-element array indicating the ordering of the Longitude, Latitude, and radial dimensions. The array must contain some permutation of the set (0,1,2). The permuation (0,1,2) indicates that longitude is the fastest varying dimesion, then latitude then radius. The permutation (2,1,0) indicates that radius is fastest, then latitude, and so on.
[in]periodicA three-element boolean vector indicating which i,j,k indecies, respectively, are periodic.
[in]blksAn array of blocks containing the sampled function. The dimensions of each block is given by bs. The number of blocks is given by the product of the terms:
(max[i]/bs[i] - min[i]/bs[i] + 1)

over i = 0..2.

Parameters
[in]coordsAn array of blocks with dimension and number the same as blks specifying the varying dimension grid point coordinates.
[in]varying_dimAn enumerant indicating which axis is the varying dimension: 0 for I, 1 for J, 2 for K

◆ SphericalGrid() [2/2]

VAPoR::SphericalGrid::SphericalGrid ( const size_t  bs[3],
const size_t  min[3],
const size_t  max[3],
const double  extents[6],
const size_t  permutation[3],
const bool  periodic[3],
float **  blks,
float  missing_value 
)

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

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

See also
GetValue()

Member Function Documentation

◆ CartToSph()

static void VAPoR::SphericalGrid::CartToSph ( double  x,
double  y,
double  z,
double *  phi,
double *  theta,
double *  r 
)
inlinestatic

◆ GetBoundingBox()

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

◆ GetIJKIndex()

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

◆ GetIJKIndexFloor()

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

◆ GetUserCoordinates()

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

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

◆ GetUserExtents()

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

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

Return extents in Cartesian coordinates

Definition at line 92 of file SphericalGrid.h.

◆ GetValue()

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

Get the reconstructed value of the sampled scalar function

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

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

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

Parameters
[in]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 from VAPoR::Grid.

◆ InsideGrid()

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

Return true if the specified point lies inside the grid

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

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

◆ SphToCart()

static void VAPoR::SphericalGrid::SphToCart ( double  phi,
double  theta,
double  r,
double *  x,
double *  y,
double *  z 
)
inlinestatic

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