Toast++
1.0.2 (r.539)
Forward and inverse modelling in optical tomography
|
Base class for finite element types. More...
#include <element.h>
Public Member Functions | |
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. | |
virtual Element * | Copy ()=0 |
Create a copy of the element and return a pointer to it. | |
virtual void | Initialise (const NodeList &nlist) |
Element initialisation. More... | |
virtual void | PostInitialisation (const NodeList &nlist) |
Element setup after mesh initialisation. More... | |
void | operator= (const Element &el) |
Element assignment. More... | |
virtual BYTE | Type () const =0 |
Returns an element type identifier. More... | |
virtual BYTE | VtkType () const |
Returns the VTK element type identifier, or 0 if the element doesn't have a VTK representation. | |
virtual unsigned long | GetCaps () const =0 |
Returns element capability flags. More... | |
virtual int | Dimension () const =0 |
Returns the spatial dimension of the element. More... | |
virtual int | nNode () const =0 |
Returns the number of nodes associated with the element. More... | |
virtual int | nSide () const =0 |
Returns the number of element sides. More... | |
virtual int | nSideNode (int side) const =0 |
Returns the number of vertices associated with a side. More... | |
virtual int | SideNode (int side, int node) const =0 |
Returns relative node index for a side vertex. 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 Point | Local (const NodeList &nlist, const Point &glob) const =0 |
Maps a point from global to local element coordinates. More... | |
virtual Point | NodeLocal (int node) const =0 |
Returns the local coordinates of an element node. More... | |
virtual Point | SurfToLocal (int side, const Point &p) const |
Maps a point from surface coordinates to local element coordinates. 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... | |
Element * | SideNeighbour (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 RVector | DirectionCosine (int side, RDenseMatrix &jacin)=0 |
Returns the direction cosines of a side normal. More... | |
virtual const RVector & | LNormal (int side) const =0 |
Returns a side normal in local coordinates. More... | |
virtual double | Size () const =0 |
Returns the element size. More... | |
virtual double | SideSize (int side, const NodeList &nlist) const |
Returns the size of an element side. More... | |
virtual double | DetJ (const Point &loc, const NodeList *nlist=0) const |
Returns determinant of Jacobian at a given point inside the element in the local frame. More... | |
virtual bool | LContains (const Point &loc, bool pad=true) const =0 |
Checks if a local point coordinate is inside the element. More... | |
virtual bool | GContains (const Point &glob, const NodeList &nlist) const |
Checks if a global point coordinate is inside the element. 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 | LocalShapeF (const Point &loc) const =0 |
Returns the values of the shape functions at a local point. More... | |
virtual RDenseMatrix | LocalShapeD (const Point &loc) const =0 |
Returns the values of the shape function derivatives at a local point. More... | |
virtual RVector | GlobalShapeF (const NodeList &nlist, const Point &glob) const |
Returns the values of the shape functions at a global point. More... | |
virtual RDenseMatrix | GlobalShapeD (const NodeList &nlist, const Point &glob) const |
Returns the values of the shape function derivatives at a global point. More... | |
virtual int | QuadRule (int order, const double **wght, const Point **absc) const |
Returns the weights and abscissae of quadrature rules over the element. More... | |
virtual double | IntF (int i) const =0 |
Integral of a shape function over the element. More... | |
virtual RSymMatrix | IntFF () const =0 |
Integrals of all products of two shape functions over the element. More... | |
virtual double | IntFF (int i, int j) const =0 |
Integral of a product of two shape functions over the element. More... | |
virtual double | IntFFF (int i, int j, int k) const =0 |
Integral of a product of three shape functions over the element. More... | |
virtual RSymMatrix | IntPFF (const RVector &P) const =0 |
Integrals of all products of two shape functions and a nodal function over the element. More... | |
virtual double | IntPFF (int i, int j, const RVector &P) const =0 |
Integral of a product of two shape functions and a nodal function over the element. More... | |
virtual RSymMatrix | IntDD () const =0 |
Integrals of all products of two shape function derivatives over the element. More... | |
virtual double | IntDD (int i, int j) const =0 |
Integral of a product of two shape function derivatives over the element. More... | |
virtual RVector | IntFD (int i, int j) const |
Integral of a product of a shape function and a shape function derivative over the element. More... | |
virtual double | IntFDD (int i, int j, int k) const =0 |
Integral of a product of a shape function and two shape function derivatives over the element. More... | |
virtual RSymMatrix | IntPDD (const RVector &P) const =0 |
All integrals of products of a nodal function and two shape function derivatives over the element. More... | |
virtual double | IntPDD (int i, int j, const RVector &P) const =0 |
Integrals of a product of a nodal function and two shape function derivatives over the element. More... | |
virtual RVector | BndIntF () const |
Boundary integral of all shape functions over all boundary sides of the element. More... | |
virtual double | BndIntFSide (int i, int sd) |
Surface integral of a shape function over one of the sides of the element. More... | |
virtual RSymMatrix | BndIntFF () const =0 |
Boundary integral of all products of two shape functions over all boundary sides of the element. More... | |
virtual double | BndIntFF (int i, int j)=0 |
Boundary integral of a product of two shape functions over all boundary sides of the element. More... | |
virtual double | BndIntFFSide (int i, int j, int sd)=0 |
Surface integral of a product of two shape functions over one of the sides of the element. More... | |
virtual RSymMatrix | BndIntPFF (const RVector &P) const =0 |
Surface integrals of all products of a nodal function and two shape functions over all boundary sides of the element. More... | |
virtual double | BndIntPFF (int i, int j, const RVector &P) const =0 |
Surface integrals of a product of a nodal function and two 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 | 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... | |
virtual 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... | |
virtual RSymMatrix | Intdd () const |
Integral of the product of two partial shape function derivatives over the element. More... | |
virtual 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... | |
virtual double | IntPdd (const RVector &p, int j, int k, int l, int m) const |
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 |
virtual double | IntUnitSphereP (const NodeList &nlist, const RVector &P) const |
virtual double | IntUnitSphereFF (const NodeList &nlist, const int i, const int j) const |
virtual double | IntUnitSpherePFF (const NodeList &nlist, const int i, const int j, const RVector &P) const |
virtual RVector | BndIntFX (int side, double(*func)(const Point &), const NodeList &nlist) const |
virtual RVector | BndIntFCos (int side, const Surface *surf, const RVector &cntcos, double sigma, double sup, const NodeList &nlist) const |
virtual RVector | BndIntFCos (int side, const RVector &cntcos, double a, const NodeList &nlist) const |
virtual RVector | BndIntFDelta (int side, const Surface *surf, const RVector &pos, const NodeList &nlist) 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 int | GetLocalSubsampleAbsc (const Point *&absc) const |
virtual int | GetBndSubsampleAbsc (int side, const Point *&absc) 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) |
virtual int | GlobalIntersection (const NodeList &nlist, const Point &p1, const Point &p2, Point **list)=0 |
virtual int | Intersection (const Point &, const Point &, Point **)=0 |
virtual RDenseMatrix | LocaltoGlobalMat () const |
virtual RDenseMatrix | GlobaltoLocalMat () const |
virtual RDenseMatrix | FTAMat () const |
int | Region () const |
void | SetRegion (int nr) |
virtual void | SplitSide (Mesh *mesh, int side, int newnode, Element *nbr1, Element *nbr2, Element *el1, Element *el2) |
virtual void | Bisect (Mesh *mesh, int side, int newnode, Element *nbr1, Element *nbr2, Element *el1, Element *el2) |
virtual void | MergeAndResplit (Mesh *mesh, int side, int newnode, Element *nbr1, Element *nbr2, Element *el1, Element *el2) |
Public Attributes | |
int * | Node |
ElementSubdivisionData * | subdivdata |
Element ** | sdnbhr |
int * | sdnbhridx |
Protected Attributes | |
bool * | bndside |
bool | bndel |
bool | interfaceel |
int | region |
Friends | |
class | ElementList |
class | Mesh |
std::istream & | operator>> (std::istream &i, Element &el) |
std::ostream & | operator<< (std::ostream &o, const Element &el) |
Base class for finite element types.
This class defines the common interface for all Toast FEM elements, in particular a set of integrals of combinations of shape function, shape function derivatives and nodal functions over the element (Int*).
|
inlinevirtual |
Boundary integral of all shape functions over all boundary sides of the element.
where the integration is performed over all sides of teh element that are part of the mesh surface.
Reimplemented in Tetrahedron4.
|
pure virtual |
Boundary integral of all products of two shape functions over all boundary sides of the element.
where the integration is performed over all sides of the element that are part of the mesh surface.
Implemented in Element_Unstructured, Pixel4, Voxel8, and Voxel27.
|
pure virtual |
Boundary integral of a product of two shape functions over all boundary sides of the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
where the integration is performed over all sides of the element that are part of the mesh surface.
i
and j
are boundary nodes. Implemented in Element_Unstructured, Pixel4, Voxel8, and Voxel27.
|
pure virtual |
Surface integral of a product of two shape functions over one of the sides of the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
sd | side index (range 0 .. nSide-1) |
where the integration is performed over side sd.
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Pixel4, Wedge6, Triangle6, Voxel8, Triangle10, Voxel27, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
|
inlinevirtual |
Surface integral of a shape function over one of the sides of the element.
i | node index (range 0 .. nNode-1) |
sd | side index (range 0 .. nSide-1) |
where the integration is performed over side sd.
Reimplemented in Triangle3, Tetrahedron4, Triangle6, and Triangle10.
|
pure virtual |
Surface integrals of all products of a nodal function and two shape functions over all boundary sides of the element.
P | nodal function coefficients over the mesh |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Pixel4, Wedge6, Voxel8, Triangle6, Triangle10, Voxel27, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
|
pure virtual |
Surface integrals of a product of a nodal function and two shape functions over all boundary sides of the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
P | nodal function coefficients over the mesh |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Pixel4, Wedge6, Triangle6, Voxel8, Triangle10, Voxel27, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
|
virtual |
Returns a list of boundary sides.
nlist | mesh node list |
address | of list to be filled with boundary side indices |
Returns determinant of Jacobian at a given point inside the element in the local frame.
loc | evaluation point in local coordinates |
nlist | node list of the associated mesh |
Reimplemented in Triangle3, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, and Tetrahedron4.
|
pure virtual |
Returns the spatial dimension of the element.
Implemented in Element_Unstructured_3D, Element_Unstructured_2D, Element_Structured_3D, Element_Structured_2D, Triangle3D6, and Triangle3D3.
|
pure virtual |
Returns the direction cosines of a side normal.
side | side index (>= 0) |
jacin | inverse of Jacobian |
Implemented in Triangle3, Line2D2, Triangle3D6, Tetrahedron4, Triangle10, Triangle3D3, Triangle6, Voxel8, Triangle3old, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, and Wedge6.
|
virtual |
Checks if a global point coordinate is inside the element.
glob | point in global coordinates |
nlist | mesh node list |
Reimplemented in Element_Unstructured, Triangle3, Triangle3D6, Line2D2, Triangle3D3, Triangle10, Triangle3old, Triangle6, Voxel8, Voxel27, Pixel4, Triangle10_ip, Triangle6_ip, and Quad4.
|
pure virtual |
Returns element capability flags.
Implemented in Triangle3, Tetrahedron4, Triangle10, Triangle6, Voxel8, Tetrahedron10_ip, Tetrahedron10, Voxel27, Line2D2, Pixel4, Triangle10_ip, Triangle6_ip, Quad4, Triangle3old, and Wedge6.
Maps a point from local element coordinates to global coordinates.
nlist | mesh node list |
loc | point in local element coordinates |
|
inlinevirtual |
Returns the values of the shape function derivatives at a global point.
glob | global point coordinates |
Reimplemented in Triangle3, Triangle3D6, Line2D2, Triangle3D3, Tetrahedron4, Triangle3old, Voxel8, Triangle6, Voxel27, Pixel4, and Tetrahedron10.
Referenced by Pixel4::GlobalShapeD(), and Voxel8::GlobalShapeD().
|
inlinevirtual |
Returns the values of the shape functions at a global point.
nlist | mesh node list |
glob | point coordinates in the local element system |
Reimplemented in Triangle3, Triangle3D6, Line2D2, Triangle3D3, Tetrahedron4, Triangle3old, Triangle10, Triangle6, Voxel8, Voxel27, Pixel4, and Tetrahedron10.
Referenced by Pixel4::GlobalShapeF(), and Voxel8::GlobalShapeF().
|
inline |
Checks if the element contains any surface sides.
|
inline |
Checks if the element contains any internal interface sides.
|
virtual |
Element initialisation.
Calculation of geometric element parameters, including element size and pre-calculation of geometry dependent element matrices.
nlist | Node list, containing node geometry |
Reimplemented in Element_Unstructured, Triangle3, Tetrahedron4, Triangle10, Triangle6, Triangle3D6, Tetrahedron10_ip, Voxel8, Voxel27, Tetrahedron10, Triangle10_ip, Pixel4, Triangle6_ip, Line2D2, Triangle3D3, Quad4, Triangle3old, Wedge6, and Wedge18inf.
|
inlinevirtual |
Integral of partial shape function derivative over the element.
i | node index (range 0 .. nNode-1) |
k | coordinate index for the derivative (range 0 .. Dimension-1) |
Reimplemented in Tetrahedron4.
|
pure virtual |
|
pure virtual |
Integral of a product of two shape function derivatives over the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
Implemented in Element_Unstructured, Voxel8, Pixel4, and Voxel27.
|
inlinevirtual |
Integral of the product of two partial shape function derivatives over the element.
The dimension of the returned matrix is nd x nd, where n is the number of nodes, and d is the dimension of the element (2 or 3). The matrix contains n x n blocks, where each block is of dimension d x d.
Reimplemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Triangle6, Voxel8, Pixel4, and Voxel27.
|
pure virtual |
Integral of a shape function over the element.
i | node index (range: 0 .. nNode-1) |
Implemented in Triangle3, Line2D2, Tetrahedron4, Triangle3old, Triangle10, Voxel8, Triangle6, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, and Wedge6.
|
inlinevirtual |
Integral of a product of a shape function and a shape function derivative over the element.
i | node index for shape function |
j | node index for shape function derivative |
Reimplemented in Triangle3, Line2D2, Triangle3old, and Triangle6.
|
inlinevirtual |
Integral of the product of a shape function and a partial shape function derivative over the element.
i | node index for shape function (range 0 .. nNode-1) |
j | node index for shape function derivative (range 0 .. nNode-1) |
k | coordinate index for derivative (range 0 .. Dimension-1) |
Reimplemented in Triangle3, Line2D2, Voxel8, Triangle6, Pixel4, Tetrahedron4, and Voxel27.
|
pure virtual |
Integral of a product of a shape function and two shape function derivatives over the element.
i | node index for shape function |
j | first node index for shape function derivative |
k | second node index for shape function derivative |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Voxel8, Triangle6, Triangle10, Pixel4, Voxel27, Wedge6, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
|
inlinevirtual |
Integral of the product of a shape function and two partial shape function derivatives over the element.
i | node index for shape function (range 0 .. nNode-1) |
j | node index for shape function derivative (range 0 .. nNode-1) |
k | node index for shape function derivative (range 0 .. nNode-1) |
l | coordinate index for derivative (range 0 .. Dimension-1) |
m | coordinate index for derivative (range 0 .. Dimension-1) |
Reimplemented in Triangle3, Line2D2, Voxel8, Tetrahedron4, Pixel4, and Voxel27.
|
pure virtual |
Integrals of all products of two shape functions over the element.
Implemented in Triangle3, Line2D2, Tetrahedron4, Triangle3old, Triangle10, Voxel8, Voxel27, Pixel4, Triangle6, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, Wedge6, and Wedge18inf.
|
pure virtual |
Integral of a product of two shape functions over the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Triangle10, Triangle6, Voxel8, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, and Wedge6.
|
pure virtual |
Integral of a product of three shape functions over the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
k | third node index (range 0 .. nNode-1) |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Triangle10, Voxel8, Triangle6, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, and Wedge6.
|
inlinevirtual |
Integral of the product of a nodal function and a partial shape function derivative over the element.
P | array of nodal coefficients defining the function |
j | node index for shape function derivative (range 0 .. nNode-1) |
k | coordinate index for derivative (range 0 .. Dimension-1) |
Reimplemented in Triangle3, Line2D2, Triangle6, and Tetrahedron4.
|
pure virtual |
All integrals of products of a nodal function and two shape function derivatives over the element.
P | nodal function coefficients over the mesh |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Voxel8, Pixel4, Triangle6, Wedge6, Triangle10, Voxel27, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
Referenced by Pixel4::IntPDD().
|
pure virtual |
Integrals of a product of a nodal function and two shape function derivatives over the element.
i | first node index |
j | second node index |
P | nodal function coefficients over the mesh |
Implemented in Triangle3, Line2D2, Triangle3old, Wedge6, Tetrahedron4, Voxel8, Triangle6, Triangle10, Pixel4, Voxel27, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
|
pure virtual |
Integrals of all products of two shape functions and a nodal function over the element.
P | nodal basis coefficients of a function defined over the mesh |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Triangle10, Voxel8, Triangle6, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Wedge6, and Triangle6_ip.
|
pure virtual |
Integral of a product of two shape functions and a nodal function over the element.
i | first node index (range 0 .. nNode-1) |
j | second node index (range 0 .. nNode-1) |
P | nodal basis coefficients of a function defined over the mesh |
Implemented in Triangle3, Line2D2, Triangle3old, Tetrahedron4, Triangle10, Triangle6, Voxel8, Voxel27, Pixel4, Wedge6, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, and Triangle6_ip.
|
inline |
Checks if a side is on the mesh surface.
side | side index (>= 0) |
|
virtual |
Checks if a node is part of the element.
node | absolute node index (>= 0) |
int Element::IsSide | ( | int | nn, |
int * | nd | ||
) |
Checks if a face defined by vertex nodes is part of the element.
nn | number of node indices in the supplied array |
nd | array of absolute node indices of length nn |
|
virtual |
|
pure virtual |
Checks if a local point coordinate is inside the element.
loc | point in local coordinates |
pad | Apply padding to the element to ensure that boundary points are considered inside. |
Implemented in Triangle3, Line2D2, Tetrahedron4, Triangle10, Triangle6, Triangle3old, Voxel8, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, Quad4, and Wedge6.
|
pure virtual |
Returns a side normal in local coordinates.
side | side index (>= 0) |
Implemented in Triangle3, Line2D2, Tetrahedron4, Triangle10, Triangle6, Triangle3old, Voxel8, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, and Wedge6.
Referenced by Pixel4::DirectionCosine(), and Voxel8::DirectionCosine().
Maps a point from global to local element coordinates.
nlist | mesh node list |
glob | global point coordinates |
Implemented in Triangle3, Tetrahedron4, Triangle3D6, Triangle10, Triangle6, Line2D2, Triangle3D3, Voxel8, Tetrahedron10_ip, Voxel27, Pixel4, Tetrahedron10, Triangle10_ip, Triangle3old, Triangle6_ip, Quad4, and Wedge6.
Referenced by Pixel4::Local(), and Voxel8::Local().
|
pure virtual |
Returns the values of the shape function derivatives at a local point.
loc | point coordinates in the local element system |
Implemented in Triangle3, Line2D2, Tetrahedron4, Triangle10, Triangle3old, Triangle6, Voxel8, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, Quad4, and Wedge6.
Returns the values of the shape functions at a local point.
loc | point coordinates in the local element system |
Implemented in Triangle3, Line2D2, Tetrahedron4, Triangle10, Triangle3old, Triangle6, Voxel8, Voxel27, Pixel4, Tetrahedron10, Tetrahedron10_ip, Triangle10_ip, Triangle6_ip, Quad4, and Wedge6.
|
virtual |
Translates a point to an element surface.
[in] | side | side index (>= 0) |
[in,out] | loc | point to be transformed. |
Reimplemented in Tetrahedron4.
|
pure virtual |
Returns the number of nodes associated with the element.
Implemented in Triangle3, Tetrahedron4, Triangle10, Triangle6, Voxel8, Tetrahedron10_ip, Line2D2, Tetrahedron10, Voxel27, Pixel4, Triangle10_ip, Triangle6_ip, Triangle3old, Quad4, Wedge6, and Wedge18inf.
|
pure virtual |
Returns the local coordinates of an element node.
node | node index (>= 0) |
Implemented in Triangle3, Tetrahedron4, Triangle10, Line2D2, Triangle6, Voxel8, Voxel27, Pixel4, Tetrahedron10_ip, Tetrahedron10, Triangle3old, Triangle10_ip, Triangle6_ip, Quad4, and Wedge6.
|
pure virtual |
Returns the number of element sides.
Implemented in Triangle3, Tetrahedron4, Triangle10, Triangle6, Line2D2, Voxel8, Tetrahedron10_ip, Tetrahedron10, Voxel27, Pixel4, Triangle10_ip, Triangle6_ip, Triangle3old, Quad4, and Wedge6.
|
pure virtual |
Returns the number of vertices associated with a side.
side | side index (>= 0) |
Implemented in Triangle3, Tetrahedron4, Triangle10, Triangle6, Line2D2, Voxel8, Tetrahedron10_ip, Tetrahedron10, Voxel27, Pixel4, Triangle10_ip, Triangle6_ip, Triangle3old, Quad4, and Wedge6.
void Element::operator= | ( | const Element & | el | ) |
Element assignment.
el | right-hand side element operand |
|
inlinevirtual |
Element setup after mesh initialisation.
This method is called after all mesh elements have been initialised. It can be used to perform any element initialisation steps that require the rest of the mesh to be initialised.
Reimplemented in Element_Unstructured, and Voxel8.
|
inlinevirtual |
Returns the weights and abscissae of quadrature rules over the element.
order | quadrature order (>= 1) |
wght | pointer to an array of quadrature weights |
absc | pointer to an array of quadrature abscissae |
Reimplemented in Triangle3, and Triangle3old.
|
virtual |
Returns the centre point of a side.
side | side index (>= 0) |
Element* Element::SideNeighbour | ( | int | side | ) | const |
returns a pointer to the element connected at 'side', or NULL if this is a boundary side.
side | side index (>= 0) |
int Element::SideNeighbourIndex | ( | int | side | ) | const |
Returns the list index of the element connected at 'side', or -1 if this is a boundary side.
side | side index (>= 0) |
|
pure virtual |
Returns relative node index for a side vertex.
side | side index (>= 0) |
node | side node index (>= 0) |
Implemented in Triangle3, Tetrahedron4, Triangle10, Triangle6, Line2D2, Voxel8, Tetrahedron10_ip, Tetrahedron10, Voxel27, Pixel4, Triangle10_ip, Triangle6_ip, Triangle3old, Quad4, and Wedge6.
|
inlinevirtual |
Returns the size of an element side.
side | side index (>= 0) |
nlist | mesh node list |
Reimplemented in Triangle3, Tetrahedron4, Line2D2, Tetrahedron10, and Triangle3old.
|
pure virtual |
Returns the element size.
Implemented in Element_Unstructured, Voxel8, Voxel27, and Pixel4.
int Element::SubdivisionLevel | ( | ) | const |
Returns the element's subdivison level.
Maps a point from surface coordinates to local element coordinates.
side | side index (>= 0) |
p | surface point coordinates. |
Reimplemented in Triangle3, Line2D2, Triangle3D6, Tetrahedron4, Triangle3D3, Triangle3old, and Tetrahedron10.
|
pure virtual |
Returns an element type identifier.
Implemented in Triangle3, Tetrahedron4, Triangle10, Triangle6, Triangle3D6, Voxel8, Tetrahedron10_ip, Voxel27, Tetrahedron10, Triangle10_ip, Pixel4, Triangle6_ip, Line2D2, Triangle3D3, Quad4, Triangle3old, Wedge6, and Wedge18inf.