Toast++
1.0.2 (r.539)
Forward and inverse modelling in optical tomography
|
Generic blob basis representation. More...
#include <raster_bl.h>
Public Member Functions | |
Raster_Blob (const IVector &_bdim, const IVector &_gdim, Mesh *mesh, double _sup, RDenseMatrix *bb=0) | |
Blob basis constructor. More... | |
virtual double | Value (const Point &p, int i) const |
Returns the value of a basis function at a given point. 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 | NodeValues (int node, RVector &nv) const |
Fill return all basis function values at a given mesh node. More... | |
RVector | NodeValues (int node) const |
void | Map_MeshToBasis (const RVector &mvec, RVector &bvec) const |
Map a real-valued field from mesh to basis representation. More... | |
void | Map_MeshToBasis (const CVector &mvec, CVector &bvec) const |
Map a complex-valued field from mesh to basis representation. More... | |
void | Map_MeshToGrid (const RVector &mvec, RVector &gvec) const |
Map a real-valued field from mesh to grid representation. More... | |
void | Map_MeshToGrid (const CVector &mvec, CVector &gvec) const |
Map a complex-valued field from mesh to grid representation. More... | |
void | Map_MeshToSol (const RVector &mvec, RVector &svec) const |
Map a real-valued field from mesh to solution representation. More... | |
void | Map_MeshToSol (const CVector &mvec, CVector &svec) const |
Map a complex-valued field from mesh to solution representation. More... | |
void | Map_BasisToMesh (const RVector &bvec, RVector &mvec) const |
Map a real-valued field from basis to mesh representation. More... | |
void | Map_BasisToMesh (const CVector &bvec, CVector &mvec) const |
Map a complex-valued field from basis to mesh representation. More... | |
void | Map_GridToMesh (const RVector &gvec, RVector &mvec) const |
Map a real-valued field from grid to mesh representation. More... | |
void | Map_GridToMesh (const CVector &gvec, CVector &mvec) const |
Map a complex-valued field from grid to mesh representation. More... | |
void | Map_SolToMesh (const RVector &svec, RVector &mvec) const |
Map a real-valued field from solution to mesh representation. More... | |
void | Map_SolToMesh (const CVector &svec, CVector &mvec) const |
Map a complex-valued field from solution to mesh representation. 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_SolToBasis (const RVector &svec, RVector &bvec) const |
Map a real-valued field from solution to basis representation. More... | |
void | Map_SolToBasis (const CVector &svec, CVector &bvec) const |
Map a complex-valued field from solution to basis representation. More... | |
void | Map_BasisToSol (const RVector &bvec, RVector &svec) const |
Map a real-valued field from basis to solution representation. More... | |
void | Map_BasisToSol (const CVector &bvec, CVector &svec) const |
Map a complex-valued field from basis to solution representation. More... | |
void | Map_SolToGrid (const RVector &svec, RVector &gvec) const |
Map a real-valued field from solution to grid representation. More... | |
void | Map_SolToGrid (const CVector &svec, CVector &gvec) const |
Map a complex-valued field from solution to grid representation. More... | |
void | Map_GridToSol (const RVector &gvec, RVector &svec) const |
Map a real-valued field from grid to solution representation. More... | |
void | Map_GridToSol (const CVector &gvec, CVector &svec) const |
Map a complex-valued field from grid to solution 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 IVector & | BDim () const |
Return native basis dimensions. More... | |
const IVector & | GDim () const |
Return high-res grid dimensions. More... | |
const RVector & | GSize () 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 Mesh & | mesh () 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 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 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_SolToBasis (const Solution &ssol, Solution &bsol, bool mapall=false) const |
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 Solution &msol, Solution &ssol, bool mapall=false) const |
Map a set of fields from mesh to solution representation. More... | |
const RGenericSparseMatrix & | Mesh2GridMatrix () 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 IVector & | GetGridIndices (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 |
RVector * | ImageJet (const RVector &x, double sd, bool *iflags) const |
RVector | ImageGradient (const RVector &x, double sd) const |
Protected Member Functions | |
double | ComputeBasisScale () |
void | ComputeNodeValues () |
Protected Attributes | |
RCompRowMatrix * | nodevals |
double | sup |
RVector | igrid |
int | npad |
Point | bbmin_pad |
Point | bbmax_pad |
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 | |
Mesh * | meshptr |
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... | |
RGenericSparseMatrix * | B |
transformation matrix mesh->grid | |
RGenericSparseMatrix * | BI |
transformation matrix grid->mesh | |
Generic blob basis representation.
Raster_Blob::Raster_Blob | ( | const IVector & | _bdim, |
const IVector & | _gdim, | ||
Mesh * | mesh, | ||
double | _sup, | ||
RDenseMatrix * | bb = 0 |
||
) |
Blob basis constructor.
_bdim | basis dimensions (the length of this vector is 2 or 3, corresponding to the mesh dimension). |
_gdim | High-res grid dimensions (the length of this vector is 2 or 3, corresponding to the mesh dimension). |
mesh | pointer to mesh instance |
_sup | blob support radius in units of grid spacing |
bb | bounding box for the basis representation (size 2x3 or 3x3) |
Map a real-valued field from basis to grid representation.
[in] | bvec | field in native basis representation |
[out] | gvec | field in grid representation |
Implements Raster.
Map a complex-valued field from basis to grid representation.
[in] | bvec | field in native basis representation |
[out] | gvec | field in grid representation |
Implements Raster.
Map a real-valued field from basis to mesh representation.
[in] | bvec | field vector in basis representation |
[out] | mvec | field vector in mesh representation Maps field from basis to grid representation, and then from grid to mesh representation. |
Reimplemented from Raster.
Map a complex-valued field from basis to mesh representation.
[in] | bvec | field vector in basis representation |
[out] | mvec | field vector in mesh representation Maps field from basis to grid representation, and then from grid to mesh representation. |
Reimplemented from Raster.
Map a real-valued field from basis to solution representation.
[in] | bvec | field vector in basis representation |
[out] | svec | field vector in solution representation Uses the sol2basis mask vector to eliminate basis functions outside the domain support. |
Reimplemented from Raster.
Map a complex-valued field from basis to solution representation.
[in] | bvec | field vector in basis representation |
[out] | svec | field vector in solution representation Uses the sol2basis mask vector to eliminate basis functions outside the domain support. |
Reimplemented from Raster.
Map a real-valued field from grid to basis representation.
[in] | gvec | field in grid representation |
[out] | bvec | field in native basis representation |
Implements Raster.
Map a complex-valued field from grid to basis representation.
[in] | gvec | field in grid representation |
[out] | bvec | field in native basis representation |
Implements Raster.
Map a real-valued field from grid to mesh representation.
[in] | gvec | field vector in grid representation |
[out] | mvec | field vector in mesh representation |
Reimplemented from Raster.
Map a complex-valued field from grid to mesh representation.
[in] | gvec | field vector in grid representation |
[out] | mvec | field vector in mesh representation |
Reimplemented from Raster.
Map a real-valued field from grid to solution representation.
[in] | gvec | field vector in grid representation |
[out] | svec | field vector in solution representation Calls Map_GridToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation. |
Reimplemented from Raster.
Map a complex-valued field from grid to solution representation.
[in] | gvec | field vector in grid representation |
[out] | svec | field vector in solution representation Calls Map_GridToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation. |
Reimplemented from Raster.
Map a real-valued field from mesh to basis representation.
[in] | mvec | field vector in mesh representation |
[out] | bvec | field vector in basis representation Maps field from mesh to grid representation, and then from grid to basis representation. |
Reimplemented from Raster.
Map a complex-valued field from mesh to basis representation.
[in] | mvec | field vector in mesh representation |
[out] | bvec | field vector in basis representation Maps field from mesh to grid representation, and then from grid to basis representation. |
Reimplemented from Raster.
Map a real-valued field from mesh to grid representation.
[in] | mvec | field vector in mesh representation |
[out] | bvec | field vector in grid representation |
Reimplemented from Raster.
Map a complex-valued field from mesh to grid representation.
[in] | mvec | field vector in mesh representation |
[out] | bvec | field vector in grid representation |
Reimplemented from Raster.
Map a real-valued field from mesh to solution representation.
[in] | mvec | field vector in mesh representation |
[out] | svec | field vector in solution representation Calls Map_MeshToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation. |
Reimplemented from Raster.
Map a complex-valued field from mesh to solution representation.
[in] | mvec | field vector in mesh representation |
[out] | svec | field vector in solution representation Calls Map_MeshToBasis to create an intermediate basis vector, followed by Map_BasisToSol to map to solution representation. |
Reimplemented from Raster.
Map a real-valued field from solution to basis representation.
[in] | svec | field vector in solution representation |
[out] | bvec | field vector in basis representation Uses the sol2basis mask vector to copy voxels within the support of the domain. External voxels are set to zero. |
Reimplemented from Raster.
Map a complex-valued field from solution to basis representation.
[in] | svec | field vector in solution representation |
[out] | bvec | field vector in basis representation Uses the sol2basis mask vector to copy voxels within the support of the domain. External voxels are set to zero. |
Reimplemented from Raster.
Map a real-valued field from solution to grid representation.
[in] | svec | field vector in solution representation |
[out] | gvec | field vector in grid representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToGrid to map to grid representation. |
Reimplemented from Raster.
Map a complex-valued field from solution to grid representation.
[in] | svec | field vector in solution representation |
[out] | gvec | field vector in grid representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToGrid to map to grid representation. |
Reimplemented from Raster.
Map a real-valued field from solution to mesh representation.
[in] | svec | field vector in solution representation |
[out] | mvec | field vector in mesh representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToMesh to map to mesh representation. |
Reimplemented from Raster.
Map a complex-valued field from solution to mesh representation.
[in] | svec | field vector in solution representation |
[out] | mvec | field vector in mesh representation Calls Map_SolToBasis to create an intermediate basis vector, followed by Map_BasisToMesh to map to mesh representation. |
Reimplemented from Raster.
|
virtual |
Fill return all basis function values at a given mesh node.
[in] | node | mesh node index (>= 0) |
[out] | nv | vector of basis function values |
|
virtual |
Returns the value of a basis function at a given point.
p | point coordinates |
i | linear basis index |
Reimplemented from Raster.
|
pure virtual |
Value of basis function b_i at point p This does not check for mesh support.
Implements Raster.
Implemented in Raster_HanningBlob, Raster_BesselBlob, Raster_GaussBlob, Raster_RampBlob, and Raster_SplineBlob.