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

Pixel (or voxel) basis representation. More...

#include <raster_px2.h>

Inheritance diagram for Raster_Pixel2:
Raster

Public Member Functions

 Raster_Pixel2 (const IVector &_bdim, const IVector &_gdim, Mesh *mesh, RDenseMatrix *bb=0, double _map_tol=1e-10)
 Pixel basis constructor. More...
 
 ~Raster_Pixel2 ()
 Pixel basis destructor.
 
double Value_nomask (const Point &p, int i, bool is_solidx=true) const
 Value of basis function b_i at point p This does not check for mesh support. More...
 
void Map_GridToBasis (const RVector &gvec, RVector &bvec) const
 Map a real-valued field from grid to basis representation. More...
 
void Map_GridToBasis (const CVector &gvec, CVector &bvec) const
 Map a complex-valued field from grid to basis representation. More...
 
void Map_BasisToGrid (const RVector &bvec, RVector &gvec) const
 Map a real-valued field from basis to grid representation. More...
 
void Map_BasisToGrid (const CVector &bvec, CVector &gvec) const
 Map a complex-valued field from basis to grid representation. More...
 
void Map_MeshToBasis (const RVector &mvec, RVector &bvec) const
 Map a real-valued field from mesh to basis representation. More...
 
void Map_BasisToMesh (const RVector &bvec, RVector &mvec) const
 Map a real-valued field from basis to mesh representation. More...
 
void Map_BasisToSol (const RVector &bvec, RVector &svec) const
 Map a real-valued field from basis to solution representation. More...
 
void Map_SolToBasis (const RVector &svec, RVector &bvec) const
 Map a real-valued field from solution to basis representation. More...
 
void Map_MeshToSol (const RVector &mvec, RVector &svec) const
 Map a real-valued field from mesh to solution representation. representation. More...
 
void Map_SolToMesh (const RVector &svec, RVector &mvec) const
 Map a real-valued field from solution to mesh representation. More...
 
- Public Member Functions inherited from Raster
 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 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 Solution &gsol, Solution &bsol, bool mapall=false) const
 Map a set of fields from grid to basis 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 CVector &mvec, CVector &bvec) const
 Map a complex-valued field from mesh to basis 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 CVector &bvec, CVector &svec) const
 Map a complex-valued field from basis to solution 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 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 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
 

Additional Inherited Members

- Protected Attributes inherited from Raster
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

Pixel (or voxel) basis representation.

Constructor & Destructor Documentation

Raster_Pixel2::Raster_Pixel2 ( const IVector _bdim,
const IVector _gdim,
Mesh mesh,
RDenseMatrix bb = 0,
double  _map_tol = 1e-10 
)

Pixel basis constructor.

Parameters
_bdimbasis dimensions (the length of this vector is 2 or 3, corresponding to the mesh dimension).
_gdimNot used by this basis type. Must be identical to _bdim.
meshpointer to mesh instance
bbbounding box for the basis representation (size 2x3 or 3x3)
_map_tolTolerance limit for least-squares solver in Map_MeshToBasis and Map_BasisToMesh
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

void Raster_Pixel2::Map_BasisToGrid ( const RVector bvec,
RVector gvec 
) const
virtual

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

Parameters
[in]bvecfield in native basis representation
[out]gvecfield in grid representation
Note
If grid and basis dimensions are identical, this operator just copies bvec into gvec. Otherwise it uses a linear interpolation scheme.

Implements Raster.

void Raster_Pixel2::Map_BasisToGrid ( const CVector bvec,
CVector gvec 
) const
virtual

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

Parameters
[in]bvecfield in native basis representation
[out]gvecfield in grid representation
Note
If grid and basis dimensions are identical, this operator just copies bvec into gvec. Otherwise it uses a linear interpolation scheme.

Implements Raster.

void Raster_Pixel2::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

Reimplemented from Raster.

void Raster_Pixel2::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
Note
The 'basis' representation contains the full bounding box, while the 'solution' representation has the masked voxels removed.

Reimplemented from Raster.

void Raster_Pixel2::Map_GridToBasis ( const RVector gvec,
RVector bvec 
) const
virtual

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

Parameters
[in]gvecfield in grid representation
[out]bvecfield in native basis representation
Note
If grid and basis dimensions are identical, this operator just copies gvec into bvec. Otherwise it uses a linear interpolation scheme.

Implements Raster.

void Raster_Pixel2::Map_GridToBasis ( const CVector gvec,
CVector bvec 
) const
virtual

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

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

Implements Raster.

void Raster_Pixel2::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

Reimplemented from Raster.

void Raster_Pixel2::Map_MeshToSol ( const RVector mvec,
RVector svec 
) const
virtual

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

Parameters
[in]mvecfield vector in mesh representation
[out]svecfield vector in solution representation

Reimplemented from Raster.

void Raster_Pixel2::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
Note
The 'basis' representation contains the full bounding box, while the 'solution' representation has the masked voxels removed.
The masked voxels in bvec are set to zero.

Reimplemented from Raster.

void Raster_Pixel2::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

Reimplemented from Raster.

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

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

See Also
Value

Implements Raster.


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