Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
Public Member Functions | List of all members
Triangle3D3 Class Reference
Inheritance diagram for Triangle3D3:
Triangle3 Element_Unstructured_2D Element_Unstructured Element

Public Member Functions

 Triangle3D3 (const Triangle3D3 &el)
 
void Initialise (const NodeList &nlist)
 Initialises the triangle. More...
 
BYTE Type () const
 Returns an element type identifier. More...
 
int Dimension () const
 Returns the spatial dimension of the element. More...
 
Point Local (const NodeList &nlist, const Point &glob) const
 Maps a point from global to local element coordinates. More...
 
Point SurfToLocal (int side, const Point &p) const
 Maps a point from surface coordinates to local element coordinates. More...
 
RVector DirectionCosine (int side, RDenseMatrix &jacin)
 Returns the direction cosines of a side normal. More...
 
bool GContains (const Point &glob, const NodeList &nlist) const
 Checks if a global point coordinate is inside the element. More...
 
RVector GlobalShapeF (const NodeList &nlist, const Point &glob) const
 Returns the values of the shape functions at a global point. More...
 
RDenseMatrix GlobalShapeD (const NodeList &nlist, const Point &glob) const
 Returns the values of the shape function derivatives at a global point. More...
 
int GlobalIntersection (const NodeList &nlist, const Point &p1, const Point &p2, Point **list)
 
RVector LocalShapeQF (const Point &loc) const
 
RDenseMatrix LocalShapeQD (const Point &loc) const
 
double IntUnitSphereP (const NodeList &nlist, const RVector &P) const
 
double IntUnitSpherePFF (const NodeList &nlist, const int i, const int j, const RVector &P) const
 
double IntUnitSphereFF (const NodeList &nlist, const int i, const int j) const
 
RDenseMatrix LocaltoGlobalMat () const
 
RDenseMatrix GlobaltoLocalMat () const
 
RDenseMatrix FTAMat () const
 
- Public Member Functions inherited from Triangle3
 Triangle3 ()
 Creates a new triangle, but does not assign any node indices.
 
 Triangle3 (const Triangle3 &el)
 Creates a new triangle as a copy of an existing one. More...
 
 ~Triangle3 ()
 Destroys the triangle.
 
ElementCopy ()
 Create a copy of the element and return a pointer to it.
 
BYTE VtkType () const
 Returns the VTK element type identifier, or 0 if the element doesn't have a VTK representation.
 
unsigned long GetCaps () const
 Returns element capability flags. More...
 
int nNode () const
 Returns the number of nodes associated with the element. More...
 
int nSide () const
 Returns the number of element sides. More...
 
int nSideNode (int) const
 Returns the number of vertices associated with a side. More...
 
int SideNode (int side, int node) const
 Returns relative node index for a side vertex. More...
 
double SideSize (int sd, const NodeList &nlist) const
 Returns the size of an element side. More...
 
double DetJ (const Point &loc, const NodeList *nlist=0) const
 Returns determinant of Jacobian. More...
 
Point NodeLocal (int node) const
 Returns the local coordinates of an element node. More...
 
const RVectorLNormal (int side) const
 Returns a side normal in local coordinates. More...
 
bool LContains (const Point &loc, bool pad=true) const
 Checks if a local point coordinate is inside the element. More...
 
RVector LocalShapeF (const Point &loc) const
 Returns the values of the shape functions at a local point. More...
 
RDenseMatrix LocalShapeD (const Point &loc) const
 Returns the values of the shape function derivatives at a local point. More...
 
double IntF (int i) const
 Integral of a shape function over the element. More...
 
RSymMatrix IntFF () const
 Integrals of all products of two shape functions over the element. More...
 
int QuadRule (int order, const double **wght, const Point **absc) const
 Returns the weights and abscissae of quadrature rules over the element. More...
 
double IntFF (int i, int j) const
 Integral of a product of two shape functions over the element. More...
 
double IntFFF (int i, int j, int k) const
 Integral of a product of three shape functions over the element. More...
 
RSymMatrix IntPFF (const RVector &P) const
 Integrals of all products of two shape functions and a nodal function over the element. More...
 
double IntPFF (int i, int j, const RVector &P) const
 Integral of a product of two shape functions and a nodal function over the element. More...
 
RVector IntFD (int i, int j) const
 Integral of a product of a shape function and a shape function derivative over the element. More...
 
double IntFDD (int i, int j, int k) const
 Integral of a product of a shape function and two shape function derivatives over the element. More...
 
RSymMatrix IntPDD (const RVector &P) const
 All integrals of products of a nodal function and two shape function derivatives over the element. More...
 
double IntPDD (int i, int j, const RVector &P) const
 Integrals of a product of a nodal function and two shape function derivatives over the element. More...
 
RSymMatrix BndIntPFF (const RVector &P) const
 Surface integrals of all products of a nodal function and two shape functions over all boundary sides of the element. More...
 
double BndIntPFF (int i, int j, const RVector &P) const
 Surface integrals of a product of a nodal function and two shape functions over all boundary sides of the element. More...
 
double IntFd (int i, int j, int k) const
 Integral of the product of a shape function and a partial shape function derivative over the element. More...
 
double IntPd (const RVector &P, int j, int k) const
 Integral of the product of a nodal function and a partial shape function derivative over the element. More...
 
double IntFdd (int i, int j, int k, int l, int m) const
 Integral of the product of a shape function and two partial shape function derivatives over the element. More...
 
double IntPdd (const RVector &p, int j, int k, int l, int m) const
 
RSymMatrix Intdd () const
 Integral of the product of two partial shape function derivatives over the element. More...
 
double BndIntFSide (int i, int sd)
 Surface integral of a shape function over one of the sides of the element. More...
 
double BndIntFFSide (int i, int j, int sd)
 Surface integral of a product of two shape functions over one of the sides of the element. More...
 
RVector BndIntFX (int side, double(*func)(const Point &), const NodeList &nlist) const
 
RVector BndIntFCos (int side, const Surface *surf, const RVector &cntcos, double sigma, double sup, const NodeList &nlist) const
 
RVector BndIntFCos (int side, const RVector &cntcos, double a, const NodeList &nlist) const
 
RVector BndIntFDelta (int side, const Surface *surf, const RVector &pos, const NodeList &nlist) const
 
int GetLocalSubsampleAbsc (const Point *&absc) const
 
int GetBndSubsampleAbsc (int side, const Point *&absc) const
 
int Intersection (const Point &p1, const Point &p2, Point **pi)
 
void SplitSide (Mesh *mesh, int side, int newnode, Element *nbr1, Element *nbr2, Element *el1, Element *el2)
 
void Bisect (Mesh *mesh, int side, int newnode, Element *nbr1, Element *nbr2, Element *el1, Element *el2)
 
void MergeAndResplit (Mesh *mesh, int side, int newnode, Element *nbr1, Element *nbr2, Element *el1, Element *el2)
 
- Public Member Functions inherited from Element_Unstructured_2D
 Element_Unstructured_2D (const Element_Unstructured_2D &el)
 
- Public Member Functions inherited from Element_Unstructured
 Element_Unstructured (const Element_Unstructured &el)
 
virtual void PostInitialisation (const NodeList &nlist)
 Element setup after mesh initialisation. More...
 
double Size () const
 Returns the element size. More...
 
void operator= (const Element_Unstructured &el)
 
RSymMatrix IntDD () const
 Integrals of all products of two shape function derivatives over the element. More...
 
double IntDD (int i, int j) const
 Integral of a product of two shape function derivatives over the element. More...
 
RSymMatrix BndIntFF () const
 Boundary integral of all products of two shape functions over all boundary sides of the element. More...
 
double BndIntFF (int i, int j)
 Boundary integral of a product of two shape functions over all boundary sides of the element. More...
 
- Public Member Functions inherited from Element
 Element ()
 Creates a new element with no nodes.
 
 Element (const Element &el)
 Creates a new element as a copy of 'el'.
 
virtual ~Element ()
 Destroys the element.
 
void operator= (const Element &el)
 Element assignment. More...
 
virtual bool IsNode (int node)
 Checks if a node is part of the element. More...
 
int IsSide (int nn, int *nd)
 Checks if a face defined by vertex nodes is part of the element. More...
 
virtual bool IsSideNode (int side, int node)
 Checks if a node is associated with a specific side. More...
 
bool IsBoundarySide (int side)
 Checks if a side is on the mesh surface. More...
 
bool HasBoundarySide ()
 Checks if the element contains any surface sides. More...
 
bool HasInterfaceSide ()
 Checks if the element contains any internal interface sides. More...
 
virtual void MapToSide (int side, Point &loc) const
 Translates a point to an element surface. More...
 
virtual Point SideCentre (int side) const
 Returns the centre point of a side. More...
 
ElementSideNeighbour (int side) const
 returns a pointer to the element connected at 'side', or NULL if this is a boundary side. More...
 
int SideNeighbourIndex (int side) const
 Returns the list index of the element connected at 'side', or -1 if this is a boundary side. More...
 
Point Global (const NodeList &nlist, const Point &loc) const
 Maps a point from local element coordinates to global coordinates. More...
 
virtual RDenseMatrix Elgeom (const NodeList &nlist) const
 Returns the element's global node coordinates. More...
 
virtual int BndSideList (const NodeList &nlist, int *list)
 Returns a list of boundary sides. More...
 
void InitNeighbourSupport ()
 
void InitSubdivisionSupport ()
 
int SubdivisionLevel () const
 Returns the element's subdivison level. More...
 
virtual void Subdivide (Mesh *mesh)
 
virtual RVector BndIntF () const
 Boundary integral of all shape functions over all boundary sides of the element. More...
 
virtual double Intd (int i, int k) const
 Integral of partial shape function derivative over the element. More...
 
virtual double IntFfd (int i, int j, int k, int l) const
 
virtual double IntPfd (const RVector &p, int j, int k, int l) const
 
int GetSubsampleFD (int &n, double *&wght, Point *&absc, RVector *&F, RDenseMatrix *&D, const NodeList &nlist) const
 
int GetBndSubsampleFD (int side, int &n, double *&wght, Point *&absc, RVector *&F, RDenseMatrix *&D, const NodeList &nlist) const
 
virtual RDenseMatrix StrainDisplacementMatrix (const Point &glob) const
 
virtual RDenseMatrix ElasticityStiffnessMatrix (const RDenseMatrix &D) const
 
virtual RDenseMatrix ElasticityStiffnessMatrix (double E, double nu) const
 
RDenseMatrix ElasticStrainDisplacement (const RVector &loc, const RDenseMatrix &gder) const
 
virtual RDenseMatrix IsotropicElasticityMatrix (double E, double nu) const
 
virtual RVector InitialStrainVector (double E, double nu, const RVector &e0)
 
virtual RVector ThermalExpansionVector (double E, double nu, double alpha, double dT)
 
virtual RVector DThermalExpansionVector (double E, double nu)
 
int Region () const
 
void SetRegion (int nr)
 

Additional Inherited Members

- Public Attributes inherited from Element
int * Node
 
ElementSubdivisionDatasubdivdata
 
Element ** sdnbhr
 
int * sdnbhridx
 
- Protected Member Functions inherited from Triangle3
double ComputeSize (const NodeList &nlist) const
 
RSymMatrix ComputeIntDD (const NodeList &nlist) const
 
void ComputeIntFD (const NodeList &nlist)
 
RSymMatrix ComputeBndIntFF (const NodeList &nlist) const
 
- Protected Attributes inherited from Triangle3
double * intfd_0
 
double * intfd_1
 
double a0
 
double b0
 
double c0
 
double a1
 
double b1
 
double c1
 
double a2
 
double b2
 
double c2
 
- Protected Attributes inherited from Element_Unstructured
double size
 
RSymMatrix intdd
 
RSymMatrix intbff
 
Point bbmin
 
Point bbmax
 
- Protected Attributes inherited from Element
bool * bndside
 
bool bndel
 
bool interfaceel
 
int region
 

Member Function Documentation

int Triangle3D3::Dimension ( ) const
inlinevirtual

Returns the spatial dimension of the element.

Returns
2 for 2-D elements (triangles, pixels etc.), 3 for 3-D elements (tetrahedra, voxels, etc.)

Reimplemented from Element_Unstructured_2D.

RVector Triangle3D3::DirectionCosine ( int  side,
RDenseMatrix jacin 
)
virtual

Returns the direction cosines of a side normal.

Parameters
sideside index (>= 0)
jacininverse of Jacobian
Returns
direction cosines of the normal to side in global coordinates.
See Also
LNormal

Reimplemented from Triangle3.

bool Triangle3D3::GContains ( const Point glob,
const NodeList nlist 
) const
virtual

Checks if a global point coordinate is inside the element.

Parameters
globpoint in global coordinates
nlistmesh node list
Returns
true if the point is inside the element, false otherwise.
See Also
LContains

Reimplemented from Triangle3.

RDenseMatrix Triangle3D3::GlobalShapeD ( const NodeList nlist,
const Point glob 
) const
virtual

Returns the values of the shape function derivatives at a global point.

Parameters
globglobal point coordinates
Returns
Matrix (size nNode x Dimension) of shape function derivatives $ du_i/d x_j $ at the point for all shape functions associated with the element.
See Also
LocalShapeF, GlobalShapeF, GlobalShapeD

Reimplemented from Triangle3.

RVector Triangle3D3::GlobalShapeF ( const NodeList nlist,
const Point glob 
) const
virtual

Returns the values of the shape functions at a global point.

Parameters
nlistmesh node list
globpoint coordinates in the local element system
Returns
Vector (size nNode) of shape function values at the point for all shape functions associated with the element.
See Also
LocalShapeF, LocalShapeD, GlobalShapeD

Reimplemented from Triangle3.

void Triangle3D3::Initialise ( const NodeList nlist)
virtual

Initialises the triangle.

Parameters
nlistmesh node list
Note
This method should be called one the mesh node indices have been assigned, to allow setting up element parameters.
Calculates Natural coordinates, and subsampling abscissae. Precalculates some integrals.

Reimplemented from Triangle3.

Point Triangle3D3::Local ( const NodeList nlist,
const Point glob 
) const
virtual

Maps a point from global to local element coordinates.

Parameters
nlistmesh node list
globglobal point coordinates
Returns
Point coordinates in the element's local reference system.
Note
Point glob should have dimension Dimension()
glob does not have to be located inside the element.
See Also
NodeLocal, Global

Reimplemented from Triangle3.

Point Triangle3D3::SurfToLocal ( int  side,
const Point p 
) const
virtual

Maps a point from surface coordinates to local element coordinates.

Parameters
sideside index (>= 0)
psurface point coordinates.
Returns
Point coordinates in the element's local reference system.
Note
Surface point p should have dimension Dimension()-1
See Also
Local, NodeLocal

Reimplemented from Triangle3.

BYTE Triangle3D3::Type ( ) const
inlinevirtual

Returns an element type identifier.

Returns
Element type identifier (see Element type identifiers).

Reimplemented from Triangle3.

References ELID_TRI3D3.


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