Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
Public Member Functions | Protected Attributes | List of all members
Raster Class Referenceabstract

Base class for mapping between mesh and an independent basis representation. More...

#include <raster.h>

Inheritance diagram for Raster:
Raster_Blob Raster_CubicPixel Raster_Pixel Raster_Pixel2 Raster_BesselBlob Raster_CubicSpline Raster_GaussBlob Raster_HanningBlob Raster_RampBlob Raster_SplineBlob

Public Member Functions

 Raster (const IVector &_bdim, const IVector &_gdim, Mesh *mesh, RDenseMatrix *bb=0)
 Raster constructor. More...
 
virtual ~Raster ()
 Raster destructor.
 
const int Dim () const
 Return problem dimension. More...
 
const IVectorBDim () const
 Return native basis dimensions. More...
 
const IVectorGDim () const
 Return high-res grid dimensions. More...
 
const RVectorGSize () const
 Return size of the grid bounding box. More...
 
int GLen () const
 Return the length of the grid vector representation.
 
int BLen () const
 Return the length of the basis vector representation.
 
int SLen () const
 Return the length of the solution vector representation.
 
const Meshmesh () const
 Return the mesh associated with the basis.
 
int * Elref () const
 Returns a pointer to the element reference array for all grid points. More...
 
int * BElref () const
 Returns a pointer to the element reference array for all basis points. More...
 
void BasisVoxelPositions (RDenseMatrix &pos) const
 Returns the positions of all basis points in a matrix. More...
 
void SolutionVoxelPositions (RDenseMatrix &pos) const
 Returns the positions of all solution points in a matrix. More...
 
virtual void Sample (const RVector &bvec, const IVector &grd, RVector &img) const
 
virtual double Value (const Point &p, int i) const
 Value of basis function b_i at point p Identical to Value_nomask, but returns zero if p is outside the mesh. More...
 
virtual double Value_nomask (const Point &p, int i, bool is_solidx=true) const =0
 Value of basis function b_i at point p This does not check for mesh support. More...
 
virtual void Map_MeshToGrid (const RVector &mvec, RVector &gvec) const
 Map a real-valued field from mesh to grid representation. More...
 
virtual void Map_MeshToGrid (const CVector &mvec, CVector &gvec) const
 Map a complex-valued field from mesh to grid representation. More...
 
virtual void Map_MeshToGrid (const Solution &msol, Solution &gsol, bool mapall=false) const
 Map a set of fields from mesh to grid representation. More...
 
virtual void Map_GridToMesh (const RVector &gvec, RVector &mvec) const
 Map a real-valued field from grid to mesh representation. More...
 
virtual void Map_GridToMesh (const CVector &gvec, CVector &mvec) const
 Map a complex-valued field from grid to mesh representation. More...
 
virtual void Map_GridToMesh (const Solution &gsol, Solution &msol, bool mapall=false) const
 Map a set of fields from grid to mesh representation. More...
 
virtual void Map_GridToBasis (const RVector &gvec, RVector &bvec) const =0
 Map a real-valued field from grid to basis representation. More...
 
virtual void Map_GridToBasis (const CVector &gvec, CVector &bvec) const =0
 Map a complex-valued field from grid to basis representation. More...
 
virtual void Map_GridToBasis (const Solution &gsol, Solution &bsol, bool mapall=false) const
 Map a set of fields from grid to basis representation. More...
 
virtual void Map_BasisToGrid (const RVector &bvec, RVector &gvec) const =0
 Map a real-valued field from basis to grid representation. More...
 
virtual void Map_BasisToGrid (const CVector &bvec, CVector &gvec) const =0
 Map a complex-valued field from basis to grid representation. More...
 
virtual void Map_BasisToGrid (const Solution &bsol, Solution &gsol, bool mapall=false) const
 Map a set of fields from basis to grid representation. More...
 
virtual void Map_MeshToBasis (const RVector &mvec, RVector &bvec) const
 Map a real-valued field from mesh to basis representation. More...
 
virtual void Map_MeshToBasis (const CVector &mvec, CVector &bvec) const
 Map a complex-valued field from mesh to basis representation. More...
 
virtual void Map_BasisToMesh (const RVector &bvec, RVector &mvec) const
 Map a real-valued field from basis to mesh representation. More...
 
virtual void Map_BasisToMesh (const CVector &bvec, CVector &mvec) const
 Map a complex-valued field from basis to mesh representation. More...
 
virtual void Map_BasisToSol (const RVector &bvec, RVector &svec) const
 Map a real-valued field from basis to solution representation. More...
 
virtual void Map_BasisToSol (const CVector &bvec, CVector &svec) const
 Map a complex-valued field from basis to solution representation. More...
 
virtual void Map_SolToBasis (const RVector &svec, RVector &bvec) const
 Map a real-valued field from solution to basis representation. More...
 
virtual void Map_SolToBasis (const CVector &svec, CVector &bvec) const
 Map a complex-valued field from solution to basis representation. More...
 
virtual void Map_SolToBasis (const Solution &ssol, Solution &bsol, bool mapall=false) const
 
virtual void Map_SolToGrid (const RVector &svec, RVector &gvec) const
 Map a real-valued field from solution to grid representation. More...
 
virtual void Map_SolToGrid (const CVector &svec, CVector &gvec) const
 Map a complex-valued field from solution to grid representation. More...
 
virtual void Map_GridToSol (const RVector &gvec, RVector &svec) const
 Map a real-valued field from grid to solution representation. More...
 
virtual void Map_GridToSol (const CVector &gvec, CVector &svec) const
 Map a complex-valued field from grid to solution representation. More...
 
virtual void Map_SolToMesh (const RVector &svec, RVector &mvec) const
 Map a real-valued field from solution to mesh representation. More...
 
virtual void Map_SolToMesh (const CVector &svec, CVector &mvec) const
 Map a complex-valued field from solution to mesh representation. More...
 
virtual void Map_SolToMesh (const Solution &ssol, Solution &msol, bool mapall=false) const
 
virtual void Map_ActiveSolToMesh (const RVector &asol, Solution &msol) const
 
virtual void Map_MeshToSol (const RVector &mvec, RVector &svec) const
 Map a real-valued field from mesh to solution representation. More...
 
virtual void Map_MeshToSol (const CVector &mvec, CVector &svec) const
 Map a complex-valued field from mesh to solution representation. More...
 
virtual void Map_MeshToSol (const Solution &msol, Solution &ssol, bool mapall=false) const
 Map a set of fields from mesh to solution representation. More...
 
const RGenericSparseMatrixMesh2GridMatrix () const
 Return the mesh->grid transformation matrix.
 
virtual const
RGenericSparseMatrix
Mesh2BasisMatrix () const
 Return the mesh->basis transformation matrix. More...
 
virtual const
RGenericSparseMatrix
Basis2MeshMatrix () const
 Return the basis->mesh transformation matrix. More...
 
virtual int GetSolIdx (int basisidx) const
 Map a linear basis index into a linear solution index. More...
 
virtual int GetSolIdx (const IVector &crd) const
 Map grid indices into linear solution index. More...
 
virtual int GetBasisIdx (int solidx) const
 Map a linear solution index into a linear basis index. More...
 
virtual void GetBasisIndices (int basisidx, IVector &crd) const
 Map a linear basis index into the xy or xyz indices of the corresponding voxel. More...
 
virtual IVectorGetGridIndices (int grididx) const
 Map a linear basis index into the xy or xyz indices of the corresponding voxel. More...
 
int Basis2Sol (int i) const
 Returns the corresponding solution vector index for a given basis vector index. More...
 
int Sol2Basis (int j) const
 Returns the corresponding basis vector index for a given solution vector index. More...
 
void NeighbourGraph (idxtype *&rowptr, idxtype *&colidx, int &nzero) const
 
void NeighbourGraph (ICompRowMatrix &NG) const
 
IVector NeighbourShift (const ICompRowMatrix &NG, int i, int j) const
 
RVector RNeighbourShift (const ICompRowMatrix &NG, int i, int j) const
 
RVector SmoothImage (const RVector &x, double sd) const
 
RVectorImageJet (const RVector &x, double sd, bool *iflags) const
 
RVector ImageGradient (const RVector &x, double sd) const
 

Protected Attributes

int dim
 problem dimension (2 or 3)
 
IVector bdim
 native basis dimensions
 
IVector gdim
 high-res grid dimension
 
int blen
 length of basis vector (full bb)
 
int glen
 length of grid vector (full bb)
 
int slen
 length of solution vector (excluding mask)
 
RVector bbsize
 physical extents of grid bounding box
 
Meshmeshptr
 mesh pointer
 
Point bbmin
 
Point bbmax
 grid bounding box
 
int * belref
 basis->element reference list (length blen)
 
int * gelref
 grid->element reference list (length glen)
 
IVector paddim
 dimensions padded to a power of 2
 
int padlen
 vector length of padded grid
 
RVector bsupport
 
IVector basis2sol
 basis->solution permutation vector. More...
 
IVector sol2basis
 Solution mask vector. More...
 
RGenericSparseMatrixB
 transformation matrix mesh->grid
 
RGenericSparseMatrixBI
 transformation matrix grid->mesh
 

Detailed Description

Base class for mapping between mesh and an independent basis representation.

Derived classes implement specific basis representations.

Constructor & Destructor Documentation

Raster::Raster ( const IVector _bdim,
const IVector _gdim,
Mesh mesh,
RDenseMatrix bb = 0 
)

Raster constructor.

Parameters
_bdimnative basis dimensions (the length of this vector is 2 or 3, corresponding to the mesh dimension)
_gdimhigh-res grid dimensions (the length of this vector is 2 or 3, corresponding to the mesh dimension)
meshpointer to mesh instance
bbbounding box for the basis representation (size 2x3 or 3x3)
Note
If the bounding box is not provided, it is calculated from the mesh geometry.
The mesh pointer must remain valid for the lifetime of the basis instance.

Member Function Documentation

virtual const RGenericSparseMatrix& Raster::Basis2MeshMatrix ( ) const
inlinevirtual

Return the basis->mesh transformation matrix.

Note
Not all basis types may implement the basis->mesh transformation with an explicit matrix operator. For these classes, calling this method will throw an error.

Reimplemented in Raster_Pixel.

int Raster::Basis2Sol ( int  i) const
inline

Returns the corresponding solution vector index for a given basis vector index.

Parameters
ibasis vector index (0 .. blen-1)
Returns
Solution vector index (0 .. slen-1), or -1 if basis index is not within the support of the domain.
void Raster::BasisVoxelPositions ( RDenseMatrix pos) const

Returns the positions of all basis points in a matrix.

Parameters
[out]posmatrix of basis point positions (blen x dim)
const IVector& Raster::BDim ( ) const
inline

Return native basis dimensions.

Returns
Vector of length Dim, containing the basis basis dimensions.
int* Raster::BElref ( ) const
inline

Returns a pointer to the element reference array for all basis points.

Returns
Array of length BLen() containing the element indices (>= 0) associated with each basis point. Basis points outside the domain support are indicated by value -1.
See Also
Elref, BLen
const int Raster::Dim ( ) const
inline

Return problem dimension.

Returns
2 or 3, depending on problem dimension (corresponds to mesh dimension).
int* Raster::Elref ( ) const
inline

Returns a pointer to the element reference array for all grid points.

Returns
Array of length GLen() containing the element indices (>= 0) associated with each grid point. Grid points outside the domain support are indicated by value -1.
See Also
BElref, GLen
const IVector& Raster::GDim ( ) const
inline

Return high-res grid dimensions.

Returns
Vector of length Dim, containing the grid dimensions.
virtual int Raster::GetBasisIdx ( int  solidx) const
virtual

Map a linear solution index into a linear basis index.

Parameters
solidxsolution index (>= 0)
Returns
Linear basis index (>= 0).
Note
For a basis of dimensions nx x ny x nz, the linear basis index for voxel (i,j,k) is given by idx = i + nx*j + nx*ny*k.
The returned solution index has all masked voxels removed.
virtual void Raster::GetBasisIndices ( int  basisidx,
IVector crd 
) const
virtual

Map a linear basis index into the xy or xyz indices of the corresponding voxel.

Parameters
[in]basisidxlinear basis index (>= 0)
[out]crdVector of voxel coordinates
virtual IVector& Raster::GetGridIndices ( int  grididx) const
virtual

Map a linear basis index into the xy or xyz indices of the corresponding voxel.

Parameters
basisidxlinear basis index (>= 0)
Returns
Vector of voxel coordinates.
Note
This version is not threadsafe Map a linear grid index into the xy or xyz indices of the corresponding grid voxel.
Parameters
grididxlinear grid index (>= 0)
Returns
Vector of grid voxel coordinates.
virtual int Raster::GetSolIdx ( int  basisidx) const
virtual

Map a linear basis index into a linear solution index.

Parameters
basisidxbasis index (>= 0)
Returns
Linear solution index (>= 0) or -1 if the specified voxel is not mapped into the solution.
Note
For a basis of dimensions nx x ny x nz, the linear basis index for voxel (i,j,k) is given by idx = i + nx*j + nx*ny*k.
The returned solution index has all masked voxels removed.
virtual int Raster::GetSolIdx ( const IVector crd) const
virtual

Map grid indices into linear solution index.

Parameters
crdVector of dimension 2 or 3 with grid indices in each dimension (>= 0) \ return Linear solution index (>= 0) or -1 if the specified voxel is not mapped into the solution.
const RVector& Raster::GSize ( ) const
inline

Return size of the grid bounding box.

Returns
Vector of length Dim, containing the physical size of the grid bounding box.
Note
The grid bounding box corresponds to the mesh bounding box, unless an explicit bounding box was specified in the constructor.
virtual void Raster::Map_BasisToGrid ( const RVector bvec,
RVector gvec 
) const
pure virtual

Map a real-valued field from basis to grid representation.

Parameters
[in]bvecfield in native basis representation
[out]gvecfield in grid representation

Implemented in Raster_Pixel2, Raster_Blob, Raster_Pixel, and Raster_CubicPixel.

virtual void Raster::Map_BasisToGrid ( const CVector bvec,
CVector gvec 
) const
pure virtual

Map a complex-valued field from basis to grid representation.

Parameters
[in]bvecfield in native basis representation
[out]gvecfield in grid representation

Implemented in Raster_Pixel2, Raster_Pixel, Raster_CubicPixel, and Raster_Blob.

virtual void Raster::Map_BasisToGrid ( const Solution bsol,
Solution gsol,
bool  mapall = false 
) const
virtual

Map a set of fields from basis to grid representation.

Parameters
[in]bsolfield set in basis representation
[out]gsolfield set in grid representation
[in]mapallflag mapping all/active only fields in the set
virtual void Raster::Map_BasisToMesh ( const RVector bvec,
RVector mvec 
) const
virtual

Map a real-valued field from basis to mesh representation.

Parameters
[in]bvecfield vector in basis representation
[out]mvecfield vector in mesh representation Maps field from basis to grid representation, and then from grid to mesh representation.

Reimplemented in Raster_Pixel2, Raster_Pixel, and Raster_Blob.

virtual void Raster::Map_BasisToMesh ( const CVector bvec,
CVector mvec 
) const
virtual

Map a complex-valued field from basis to mesh representation.

Parameters
[in]bvecfield vector in basis representation
[out]mvecfield vector in mesh representation Maps field from basis to grid representation, and then from grid to mesh representation.

Reimplemented in Raster_Blob.

virtual void Raster::Map_BasisToSol ( const RVector bvec,
RVector svec 
) const
virtual

Map a real-valued field from basis to solution representation.

Parameters
[in]bvecfield vector in basis representation
[out]svecfield vector in solution representation Uses the sol2basis mask vector to eliminate basis functions outside the domain support.
Note
Derived classes should either initialise the sol2basis vector, or overload this method. If neither is done, then a default mask based on the position of each basis function (but disregarding support area) is used.

Reimplemented in Raster_Pixel2, Raster_Pixel, and Raster_Blob.

virtual void Raster::Map_BasisToSol ( const CVector bvec,
CVector svec 
) const
virtual

Map a complex-valued field from basis to solution representation.

Parameters
[in]bvecfield vector in basis representation
[out]svecfield vector in solution representation Uses the sol2basis mask vector to eliminate basis functions outside the domain support.
Note
Derived classes should either initialise the sol2basis vector, or overload this method. If neither is done, then a default mask based on the position of each basis function (but disregarding support area) is used.

Reimplemented in Raster_Blob.

virtual void Raster::Map_GridToBasis ( const RVector gvec,
RVector bvec 
) const
pure virtual

Map a real-valued field from grid to basis representation.

Parameters
[in]gvecfield in grid representation
[out]bvecfield in native basis representation

Implemented in Raster_Blob, Raster_Pixel2, Raster_Pixel, and Raster_CubicPixel.

virtual void Raster::Map_GridToBasis ( const CVector gvec,
CVector bvec 
) const
pure virtual

Map a complex-valued field from grid to basis representation.

Parameters
[in]gvecfield in grid representation
[out]bvecfield in native basis representation

Implemented in Raster_Blob, Raster_Pixel2, Raster_Pixel, and Raster_CubicPixel.

virtual void Raster::Map_GridToBasis ( const Solution gsol,
Solution bsol,
bool  mapall = false 
) const
virtual

Map a set of fields from grid to basis representation.

Parameters
[in]gsolfield set in grid representation
[out]bsolfield set in basis representation
[in]mapallflag mapping all/active only fields in the set
virtual void Raster::Map_GridToMesh ( const RVector gvec,
RVector mvec 
) const
virtual

Map a real-valued field from grid to mesh representation.

Parameters
[in]gvecfield vector in grid representation
[out]mvecfield vector in mesh representation

Reimplemented in Raster_Blob.

virtual void Raster::Map_GridToMesh ( const CVector gvec,
CVector mvec 
) const
virtual

Map a complex-valued field from grid to mesh representation.

Parameters
[in]gvecfield vector in grid representation
[out]mvecfield vector in mesh representation

Reimplemented in Raster_Blob.

virtual void Raster::Map_GridToMesh ( const Solution gsol,
Solution msol,
bool  mapall = false 
) const
virtual

Map a set of fields from grid to mesh representation.

Parameters
[in]gsolfield set in grid representation
[out]msolfield set in mesh representation
[in]mapallflag mapping all/active only fields in the set
virtual void Raster::Map_GridToSol ( const RVector gvec,
RVector svec 
) const
virtual

Map a real-valued field from grid to solution representation.

Parameters
[in]gvecfield vector in grid representation
[out]svecfield vector in solution representation Calls Map_GridToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Blob.

virtual void Raster::Map_GridToSol ( const CVector gvec,
CVector svec 
) const
virtual

Map a complex-valued field from grid to solution representation.

Parameters
[in]gvecfield vector in grid representation
[out]svecfield vector in solution representation Calls Map_GridToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Blob.

virtual void Raster::Map_MeshToBasis ( const RVector mvec,
RVector bvec 
) const
virtual

Map a real-valued field from mesh to basis representation.

Parameters
[in]mvecfield vector in mesh representation
[out]bvecfield vector in basis representation Maps field from mesh to grid representation, and then from grid to basis representation.

Reimplemented in Raster_Pixel2, Raster_Pixel, and Raster_Blob.

virtual void Raster::Map_MeshToBasis ( const CVector mvec,
CVector bvec 
) const
virtual

Map a complex-valued field from mesh to basis representation.

Parameters
[in]mvecfield vector in mesh representation
[out]bvecfield vector in basis representation Maps field from mesh to grid representation, and then from grid to basis representation.

Reimplemented in Raster_Blob.

virtual void Raster::Map_MeshToGrid ( const RVector mvec,
RVector gvec 
) const
virtual

Map a real-valued field from mesh to grid representation.

Parameters
[in]mvecfield vector in mesh representation
[out]bvecfield vector in grid representation

Reimplemented in Raster_Blob.

virtual void Raster::Map_MeshToGrid ( const CVector mvec,
CVector gvec 
) const
virtual

Map a complex-valued field from mesh to grid representation.

Parameters
[in]mvecfield vector in mesh representation
[out]bvecfield vector in grid representation

Reimplemented in Raster_Blob.

virtual void Raster::Map_MeshToGrid ( const Solution msol,
Solution gsol,
bool  mapall = false 
) const
virtual

Map a set of fields from mesh to grid representation.

Parameters
[in]msolfield set in mesh representation
[out]gsolfield set in grid representation
[in]mapallflag mapping all/active only fields in the set
virtual void Raster::Map_MeshToSol ( const RVector mvec,
RVector svec 
) const
virtual

Map a real-valued field from mesh to solution representation.

Parameters
[in]mvecfield vector in mesh representation
[out]svecfield vector in solution representation Calls Map_MeshToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Pixel2, Raster_Pixel, and Raster_Blob.

virtual void Raster::Map_MeshToSol ( const CVector mvec,
CVector svec 
) const
virtual

Map a complex-valued field from mesh to solution representation.

Parameters
[in]mvecfield vector in mesh representation
[out]svecfield vector in solution representation Calls Map_MeshToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Blob.

virtual void Raster::Map_MeshToSol ( const Solution msol,
Solution ssol,
bool  mapall = false 
) const
virtual

Map a set of fields from mesh to solution representation.

Parameters
[in]mvecfield set in mesh representation
[out]svecfield set in solution representation
[in]mapallflag mapping all/active only fields in the set Calls Map_MeshToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation.
Note
Derived classes may overload this method for a direct mapping.
virtual void Raster::Map_SolToBasis ( const RVector svec,
RVector bvec 
) const
virtual

Map a real-valued field from solution to basis representation.

Parameters
[in]svecfield vector in solution representation
[out]bvecfield vector in basis representation Uses the sol2basis mask vector to copy voxels within the support of the domain. External voxels are set to zero.
Note
Derived classes should either initialise the sol2basis vector, or overload this method. If neither is done, then a default mask based on the position of each basis function (but disregarding support area) is used.

Reimplemented in Raster_Pixel2, Raster_Pixel, and Raster_Blob.

virtual void Raster::Map_SolToBasis ( const CVector svec,
CVector bvec 
) const
virtual

Map a complex-valued field from solution to basis representation.

Parameters
[in]svecfield vector in solution representation
[out]bvecfield vector in basis representation Uses the sol2basis mask vector to copy voxels within the support of the domain. External voxels are set to zero.
Note
Derived classes should either initialise the sol2basis vector, or overload this method. If neither is done, then a default mask based on the position of each basis function (but disregarding support area) is used.

Reimplemented in Raster_Blob.

virtual void Raster::Map_SolToGrid ( const RVector svec,
RVector gvec 
) const
virtual

Map a real-valued field from solution to grid representation.

Parameters
[in]svecfield vector in solution representation
[out]gvecfield vector in grid representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToGrid to map to grid representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Blob.

virtual void Raster::Map_SolToGrid ( const CVector svec,
CVector gvec 
) const
virtual

Map a complex-valued field from solution to grid representation.

Parameters
[in]svecfield vector in solution representation
[out]gvecfield vector in grid representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToGrid to map to grid representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Blob.

virtual void Raster::Map_SolToMesh ( const RVector svec,
RVector mvec 
) const
virtual

Map a real-valued field from solution to mesh representation.

Parameters
[in]svecfield vector in solution representation
[out]mvecfield vector in mesh representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToMesh to map to mesh representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Pixel2, Raster_Pixel, and Raster_Blob.

virtual void Raster::Map_SolToMesh ( const CVector svec,
CVector mvec 
) const
virtual

Map a complex-valued field from solution to mesh representation.

Parameters
[in]svecfield vector in solution representation
[out]mvecfield vector in mesh representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToMesh to map to mesh representation.
Note
Derived classes may overload this method for a direct mapping.

Reimplemented in Raster_Blob.

virtual const RGenericSparseMatrix& Raster::Mesh2BasisMatrix ( ) const
inlinevirtual

Return the mesh->basis transformation matrix.

Note
Not all basis types may implement the basis->mesh transformation with an explicit matrix operator. For these classes, calling this method will throw an error.

Reimplemented in Raster_Pixel.

int Raster::Sol2Basis ( int  j) const
inline

Returns the corresponding basis vector index for a given solution vector index.

Parameters
jsolution vector index (0 .. slen-1)
Returns
Basis vector index (0 .. blen-1)
void Raster::SolutionVoxelPositions ( RDenseMatrix pos) const

Returns the positions of all solution points in a matrix.

Parameters
[out]posmatrix of solution point positions (slen x dim)
virtual double Raster::Value ( const Point p,
int  i 
) const
virtual

Value of basis function b_i at point p Identical to Value_nomask, but returns zero if p is outside the mesh.

See Also
Value_nomask

Reimplemented in Raster_Blob.

virtual double Raster::Value_nomask ( const Point p,
int  i,
bool  is_solidx = true 
) const
pure virtual

Value of basis function b_i at point p This does not check for mesh support.

See Also
Value

Implemented in Raster_Pixel2, Raster_Blob, Raster_Pixel, Raster_CubicPixel, Raster_HanningBlob, Raster_BesselBlob, Raster_GaussBlob, Raster_RampBlob, and Raster_SplineBlob.

Member Data Documentation

IVector Raster::basis2sol
protected

basis->solution permutation vector.

This vector is of length blen (basis dimension). Each entry contains the index of the corresponding element in the solution vector (range 0 .. slen-1), or -1 if the element is outside the support of the domain.

Note
By default, this vector is constructed from the element support array for the basis voxels. This may need to be redefined by derived classes to accommodate their basis function support areas.
RVector Raster::bsupport
protected

Vector of length blen with basis support weights. Range: 0 (no mesh support) to 1 (full mesh support) Currently only binary states (0 or 1) are supported.

IVector Raster::sol2basis
protected

Solution mask vector.

This vector is of length slen (solution dimension). Each entry contains the index of the corresponding element in the basis vector (range 0 .. blen-1). Used to map between basis representation (full bounding box) and solution representation (supported voxels only).

Note
By default, this vector is constructed from the element support array for the basis voxels. This may need to be redefined by derived classes to accommodate their basis function support areas.

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