Toast++
1.0.2 (r.539)
Forward and inverse modelling in optical tomography
|
Finite-element mesh management. More...
#include <mesh.h>
Public Member Functions | |
Mesh () | |
Constructs an empty mesh. More... | |
void | Setup (bool mark_boundary=true) |
void | Copy (const Mesh &mesh) |
int | Dimension () const |
Mesh dimension. More... | |
int | nlen () const |
Number of mesh nodes. More... | |
int | elen () const |
Number of mesh elements. More... | |
int | ilen () const |
int | nbnd () const |
RDenseMatrix | ElGeom (int el) const |
double | ElSize (int el) const |
Point | ElCentre (int el) const |
Point | ElSideCentre (int el, int sd) const |
double | ElSideSize (int el, int sd) const |
RVector | ElDirectionCosine (int el, int sd, Point *p=0) const |
double | FullSize () const |
Returns the mesh area or volume. More... | |
int | ElFind (const Point &pt) const |
int | NdFind (const Point &pt, double &dist) const |
void | Reorder (const IVector &perm) |
Re-order mesh nodes. More... | |
double | ElDist (int el1, int el2) const |
virtual void | ScaleMesh (double scale) |
virtual void | ScaleMesh (const RVector &scale) |
double | ParamAverage (const ParameterType prmtp) const |
int | MaxNodeDiff (int mode=BW_AUTO) const |
Point | NeighbourBarycentre (int node) |
void | SparseRowStructure (idxtype *&rowptr, idxtype *&colidx, int &nzero) const |
void | NeighbourCount (int *plist, int nnode, bool include_self=false) const |
void | SysMatrixStructure (int *nz, int **row_ind, int **col_ind) |
int | FindBoundarySegment (const Point &p, int *n1, int *n2, double *dist1, double *dist2) const |
double | BoundaryDistance (int node) const |
int | BoundaryList (int **bndellist, int **bndsdlist) const |
bool | PullToBoundary (const Point &p, Point &pshift, int &element, int &side, double sub=0) const |
bool | PullToBoundary_old (const Point &p, Point &pshift, int &element, int &side, double sub=0) const |
void | SetupNeighbourList () const |
void | NodeNeighbourList (int **_nnbhrs, int ***_nbhrs) const |
bool | ElConnected (int el1, int el2, int *sd1=NULL, int *sd2=NULL) |
Checks if two mesh elements have a common side. More... | |
double | Size (Point *centre=0) const |
void | BoundingBox (Point &mmin, Point &mmax, double pad=0.0) const |
Point | BndIntersect (const Point &pt1, const Point &pt2) |
void | ResetCoeff_homog (ParameterType prmtp, double val) |
void | ResetCoeff_region (ParameterType prmtp, double val, int region) |
void | ResetCoeff_sqrt (ParameterType prmtp, double cnt, double bnd) |
int | CheckConsistency () const |
void | MarkBoundary () |
RCompRowMatrix * | MassMatrix () const |
Returns the mass matrix for the mesh. More... | |
void | put (std::ostream &os, ParameterType p1=PRM_MUA, ParameterType p2=PRM_KAPPA, ParameterType p3=PRM_N) |
void | WriteVtk (ostream &os, const RVector &nim) |
Surface * | Boundary () const |
void | SetBoundary (const Surface &_boundary) |
void | PopulateNeighbourLists () |
void | InitSubdivisionSupport () |
int | RefineElement (int el) |
Refine an element by subdivision. More... | |
Public Attributes | |
NodeList | nlist |
ElementList | elist |
ParameterList | plist |
int ** | nbhrs |
int | BndType |
bool | trap_load_error |
RVector * | bnd_param |
int * | IndexBnd2Node |
int * | IndexNode2Bnd |
Friends | |
FELIB void | AddToElMatrix (const Mesh &mesh, int el, RGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToElMatrix (const Mesh &mesh, int el, FGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToElMatrix (const Mesh &mesh, int el, SCGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToElMatrix (const Mesh &mesh, int el, CGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToSysMatrix (const Mesh &mesh, RGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToSysMatrix (const Mesh &mesh, FGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToSysMatrix (const Mesh &mesh, SCGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToSysMatrix (const Mesh &mesh, CGenericSparseMatrix &M, const RVector *coeff, int mode) |
FELIB void | AddToSysMatrix (const Mesh &mesh, CGenericSparseMatrix &M, const double coeff, int mode) |
FELIB void | AddToSysMatrix (const Mesh &mesh, SCGenericSparseMatrix &M, const double coeff, int mode) |
FELIB void | AddToSysMatrix_elasticity (const Mesh &mesh, RGenericSparseMatrix &M, const RVector &modulus, const RVector &pratio) |
FELIB void | AddToRHS_elasticity (const Mesh &mesh, RVector &rhs, const RVector *coeff, int mode) |
FELIB void | AddToRHS_thermal_expansion (const Mesh &mesh, RVector &rhs, const RVector &modulus, const RVector &pratio, const RVector &thermal_expansion, double deltaT) |
FELIB void | AddElasticStrainDisplacementToSysMatrix (const Mesh &mesh, RGenericSparseMatrix &M, double E, double nu, int *nf) |
FELIB RGenericSparseMatrix * | GridMapMatrix (const Mesh &mesh, const IVector &gdim, const Point *bbmin, const Point *bbmax, const int *elref) |
Nodal to grid basis mapping. More... | |
FELIB RGenericSparseMatrix * | NodeMapMatrix (const Mesh &mesh, const IVector &gdim, const Point *bbmin, const Point *bbmax, const int *elref) |
Grid to nodal basis mapping. More... | |
FELIB RGenericSparseMatrix * | NodeMapMatrix2 (const Mesh &mesh, const IVector &gdim, const Point *bbmin, const Point *bbmax, const int *elref) |
FELIB RGenericSparseMatrix * | Grid2LinPixMatrix (const IVector &gdim, const IVector &bdim, const int *elref) |
FELIB RGenericSparseMatrix * | LinPix2GridMatrix (const IVector &gdim, const IVector &bdim, const int *elref) |
FELIB RGenericSparseMatrix * | CubicPix2GridMatrix (const IVector &bdim, const IVector &gdim, const int *elref) |
FELIB int * | GenerateElementPixelRef (const Mesh &mesh, const IVector &gdim, const Point *bbmin, const Point *bbmax) |
FELIB void | GenerateVoxelPositions (const Mesh &mesh, const IVector &gdim, const Point *bbmin, const Point *bbmax, RDenseMatrix &pos) |
FELIB void | SubsampleLinPixel (const RVector &gf, RVector &bf, const IVector &gdim, const IVector &bdim, const int *elref) |
FELIB void | RasterLinPixel (RVector &gf, const RVector &bf, const IVector &gdim, const IVector &bdim, const int *elref) |
FELIB void | SubsampleCubPixel (const RVector &gf, RVector &bf, const IVector &gdim, const IVector &bdim, const int *elref) |
FELIB void | RasterCubPixel (RVector &gf, const RVector &bf, const IVector &gdim, const IVector &bdim, const int *elref) |
FELIB std::istream & | operator>> (std::istream &i, Mesh &mesh) |
FELIB std::ostream & | operator<< (std::ostream &o, Mesh &mesh) |
Finite-element mesh management.
Defines node coordinates (via NodeList) and element connectivity (via ElementList) and provides general mesh-related methods. This version also defines node-based parameter coefficient sets for optical tomography (via ParameterList), but this should be removed.
Mesh::Mesh | ( | ) |
Constructs an empty mesh.
|
inline |
Mesh dimension.
bool Mesh::ElConnected | ( | int | el1, |
int | el2, | ||
int * | sd1 = NULL , |
||
int * | sd2 = NULL |
||
) |
Checks if two mesh elements have a common side.
[in] | el1 | element index 1 (>= 0) |
[in] | el2 | element index 2 (>= 0) |
[out] | sd1 | pointer to variable receiving the side index of the connected side in the first element |
[out] | sd2 | pointer to variable receiving the side index of the connected side in the second element |
|
inline |
Number of mesh elements.
double Mesh::FullSize | ( | ) | const |
Returns the mesh area or volume.
RCompRowMatrix* Mesh::MassMatrix | ( | ) | const |
Returns the mass matrix for the mesh.
Allocates a sparse matrix and fills it by executing the Element::IntFF operation (integral of product of two shape functions) over the elements.
|
inline |
Number of mesh nodes.
int Mesh::RefineElement | ( | int | el | ) |
Refine an element by subdivision.
el | element index (>= 0) |
void Mesh::Reorder | ( | const IVector & | perm | ) |
Re-order mesh nodes.
This method reorders the mesh nodes such that N[i] <- N[perm[i]]. The node list is rearranged, as well as the indices in the element list.
perm | Permutation list (must be of length nlen() and contain all indices from 0 to nlen()-1) |
|
friend |
Nodal to grid basis mapping.
Generate a sparse matrix which will map a nodal vector from rectangular area [bbmin,bbmax] to a regular grid of dimensions gdim.
mesh | Mesh defining the unstructured mesh basis |
gdim | Grid dimensions (# voxels) |
bbmin | min. coordinates of grid bounding box |
bbmax | max. coordinates of grid bounding box |
elref | element reference list as provided by GenerateElementPixelRef(). |
|
friend |
Grid to nodal basis mapping.
Generate a sparse matrix which will map a bitmap to a nodal image gdim contains the bitmap size, [bbmin,bbmax] is the physical area covered by the bitmap.
mesh | Mesh defining the unstructured mesh basis |
gdim | Grid dimensions (# voxels) |
bbmin | min. coordinates of grid bounding box |
bbmax | max. coordinates of grid bounding box |
elref | element reference list as provided by GenerateElementPixelRef(). |