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

Public Member Functions

 Triangle3old (const Triangle3old &el)
 
ElementCopy ()
 Create a copy of the element and return a pointer to it.
 
void Initialise (const NodeList &nlist)
 Element initialisation. More...
 
BYTE Type () const
 Returns an element type identifier. More...
 
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...
 
Point Local (const NodeList &nlist, const Point &glob) const
 Maps a point from global to local element coordinates. More...
 
Point NodeLocal (int node) const
 Returns the local coordinates of an element node. 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...
 
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...
 
bool GContains (const Point &glob, const NodeList &nlist) const
 Checks if a global 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...
 
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...
 
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...
 
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...
 
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...
 
RSymMatrix Intdd () const
 Integral of the product of two partial shape function derivatives over 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 GlobalIntersection (const NodeList &nlist, const Point &p1, const Point &p2, Point **list)
 
int Intersection (const Point &p1, const Point &p2, Point **pi)
 
- Public Member Functions inherited from Element_Unstructured_2D
 Element_Unstructured_2D (const Element_Unstructured_2D &el)
 
virtual int Dimension (void) const
 Returns the spatial dimension of the element. More...
 
- 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 BYTE VtkType () const
 Returns the VTK element type identifier, or 0 if the element doesn't have a VTK representation.
 
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 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 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 BndIntFSide (int i, int sd)
 Surface integral of a shape function over one of the 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 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
 
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)
 
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)
 

Protected Member Functions

double ComputeSize (const NodeList &nlist) const
 
RSymMatrix ComputeIntDD (const NodeList &nlist) const
 
void ComputeIntFD (const NodeList &nlist)
 
RSymMatrix ComputeBndIntFF (const NodeList &nlist) const
 

Additional Inherited Members

- Public Attributes inherited from Element
int * Node
 
ElementSubdivisionDatasubdivdata
 
Element ** sdnbhr
 
int * sdnbhridx
 
- 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

double Triangle3old::BndIntFFSide ( int  i,
int  j,
int  sd 
)
virtual

Surface integral of a product of two shape functions over one of the sides of the element.

Parameters
ifirst node index (range 0 .. nNode-1)
jsecond node index (range 0 .. nNode-1)
sdside index (range 0 .. nSide-1)
Returns
Value of the integral

\[ \int_{\partial\Omega} u_i(\vec{r}) u_j(\vec{r}) d\vec{r} \]

where the integration is performed over side sd.
Note
Nodes i and j must both belong to side sd.
See Also
BndIntFF()const, BndIntFF(int,int)

Implements Element.

RSymMatrix Triangle3old::BndIntPFF ( const RVector P) const
inlinevirtual

Surface integrals of all products of a nodal function and two shape functions over all boundary sides of the element.

Parameters
Pnodal function coefficients over the mesh
Returns
Matrix of size nNode x nNode, containing the integrals

\[ \int_{\partial\Omega} p(\vec{r}) u_i(\vec{r}) u_j(\vec{r}) d\vec{r} = \sum_k p_k \int_{\partial\Omega} u_i(\vec{r}) u_j(\vec{r}) u_k(\vec{r}) d\vec{r} \]

Note
If the element does not contain boundary sides, the returned matrix is zero.
See Also
BndIntPFF(int,int,const RVector&)const, BndIntFF

Implements Element.

double Triangle3old::BndIntPFF ( int  i,
int  j,
const RVector P 
) const
inlinevirtual

Surface integrals of a product of a nodal function and two shape functions over all boundary sides of the element.

Parameters
ifirst node index (range 0 .. nNode-1)
jsecond node index (range 0 .. nNode-1)
Pnodal function coefficients over the mesh
Returns
Value of the integral

\[ \int_{\partial\Omega} p(\vec{r}) u_i(\vec{r}) u_j(\vec{r}) d\vec{r} = \sum_k p_k \int_{\partial\Omega} u_i(\vec{r}) u_j(\vec{r}) u_k(\vec{r}) d\vec{r} \]

Note
If the element does not contain boundary sides, the returned matrix is zero.
See Also
BndIntPFF(const RVector&)const, BndIntFF

Implements Element.

RVector Triangle3old::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

Implements Element.

bool Triangle3old::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 Element_Unstructured.

unsigned long Triangle3old::GetCaps ( ) const
inlinevirtual

Returns element capability flags.

Returns
Bitflags for element caps. (see Element capability flags)

Implements Element.

References ELCAPS_SUBSAMPLING.

RDenseMatrix Triangle3old::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 Element.

RVector Triangle3old::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 Element.

void Triangle3old::Initialise ( const NodeList nlist)
virtual

Element initialisation.

Calculation of geometric element parameters, including element size and pre-calculation of geometry dependent element matrices.

Parameters
nlistNode list, containing node geometry
Note
Should be called after the node list has been generated, and whenever the node list changes.

Reimplemented from Element_Unstructured.

RSymMatrix Triangle3old::Intdd ( ) const
virtual

Integral of the product of two partial shape function derivatives over the element.

Returns
Matrix of integrals for all element nodes in each dimension.

\[ \int_\Omega \frac{\partial u_i(\vec{r})}{\partial x_j} \frac{partial u_k(\vec{r})}{\partial x_l} d\vec{r} \]

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 from Element.

double Triangle3old::IntF ( int  i) const
virtual

Integral of a shape function over the element.

Parameters
inode index (range: 0 .. nNode-1)
Returns
value of the integral

\[ \int_\Omega u_i(\vec{r}) d\vec{r} \]

See Also
IntFF, IntFFF, IntDD, IntFD

Implements Element.

RVector Triangle3old::IntFD ( int  i,
int  j 
) const
virtual

Integral of a product of a shape function and a shape function derivative over the element.

Parameters
inode index for shape function
jnode index for shape function derivative
Returns
Vector of size Dimension, containing the integral

\[ \int_\Omega u_i(\vec{r}) \nabla u_j(\vec{r}) d\vec{r} \]

See Also
IntDD, IntFDD, IntPDD

Reimplemented from Element.

double Triangle3old::IntFDD ( int  i,
int  j,
int  k 
) const
inlinevirtual

Integral of a product of a shape function and two shape function derivatives over the element.

Parameters
inode index for shape function
jfirst node index for shape function derivative
ksecond node index for shape function derivative
Returns
Value of the integral

\[ \int_\Omega u_i(\vec{r}) \nabla u_j(\vec{r}) \nabla u_k(\vec{r}) d\vec{r} \]

See Also
IntDD, IntFD, IntPDD

Implements Element.

RSymMatrix Triangle3old::IntFF ( ) const
virtual

Integrals of all products of two shape functions over the element.

Returns
Matrix of size nNode x nNode, containing the integrals

\[ \int_\Omega u_i(\vec{r}) u_j(\vec{r}) d\vec{r} \]

See Also
IntFF(int,int)const, IntF, IntFFF, IntDD, IntFD

Implements Element.

double Triangle3old::IntFF ( int  i,
int  j 
) const
virtual

Integral of a product of two shape functions over the element.

Parameters
ifirst node index (range 0 .. nNode-1)
jsecond node index (range 0 .. nNode-1)
Returns
Value of the integral

\[ \int_\Omega u_i(\vec{r}) u_j(\vec{r}) d\vec{r} \]

See Also
IntFF()const, IntF, IntFFF, IntDD, IntFD

Implements Element.

double Triangle3old::IntFFF ( int  i,
int  j,
int  k 
) const
virtual

Integral of a product of three shape functions over the element.

Parameters
ifirst node index (range 0 .. nNode-1)
jsecond node index (range 0 .. nNode-1)
kthird node index (range 0 .. nNode-1)
Returns
Value of the integral

\[ \int_\Omega u_i(\vec{r}) u_j(\vec{r}) u_k(\vec{r}) d\vec{r} \]

See Also
IntF, IntFF, IntDD, IntFD

Implements Element.

RSymMatrix Triangle3old::IntPDD ( const RVector P) const
inlinevirtual

All integrals of products of a nodal function and two shape function derivatives over the element.

Parameters
Pnodal function coefficients over the mesh
Returns
Matrix of size nNode x nNode, containing the integrals

\[ \int_\Omega p(\vec{r}) \nabla u_i(\vec{r}) \nabla u_j(\vec{r}) d\vec{r} = \sum_k p_k \int_\Omega u_k(\vec{r}) \nabla u_i(\vec{r}) \nabla u_j(\vec{r}) d\vec{r} \]

See Also
IntPDD(int,int,const RVector&)const, IntFDD, IntFD, IntPFF

Implements Element.

double Triangle3old::IntPDD ( int  i,
int  j,
const RVector P 
) const
inlinevirtual

Integrals of a product of a nodal function and two shape function derivatives over the element.

Parameters
ifirst node index
jsecond node index
Pnodal function coefficients over the mesh
Returns
Value of the integral

\[ \int_\Omega p(\vec{r}) \nabla u_i(\vec{r}) \nabla u_j(\vec{r}) d\vec{r} = \sum_k p_k \int_\Omega u_k(\vec{r}) \nabla u_i(\vec{r}) \nabla u_j(\vec{r}) d\vec{r} \]

See Also
IntPDD(const RVector&)const, IntFDD, IntFD, IntPFF

Implements Element.

RSymMatrix Triangle3old::IntPFF ( const RVector P) const
virtual

Integrals of all products of two shape functions and a nodal function over the element.

Parameters
Pnodal basis coefficients of a function defined over the mesh
Returns
Matrix of size nNode x nNode, containing the integrals

\[ \int_\Omega p(\vec{r}) u_i(\vec{r}) u_j(\vec{r}) d\vec{r} = \sum_k p_k \int_\Omega u_i(\vec{r}) u_j(\vec{r}) u_k(\vec{r}) d\vec{r} \]

See Also
IntPFF(int,int,const RVector&)const,
IntF, IntFF, IntFFF, IntFD, IntDD, IntPDD

Implements Element.

double Triangle3old::IntPFF ( int  i,
int  j,
const RVector P 
) const
virtual

Integral of a product of two shape functions and a nodal function over the element.

Parameters
ifirst node index (range 0 .. nNode-1)
jsecond node index (range 0 .. nNode-1)
Pnodal basis coefficients of a function defined over the mesh
Returns
The value of the integral

\[ \int_\Omega p(\vec{r}) u_i(\vec{r}) u_j(\vec{r}) d\vec{r} = \sum_k p_k \int_\Omega u_i(\vec{r}) u_j(\vec{r}) u_k(\vec{r}) d\vec{r} \]

See Also
IntPFF(const RVector&)const,
IntF, IntFF, IntFFF, IntFD, IntDD, IntPDD

Implements Element.

bool Triangle3old::LContains ( const Point loc,
bool  pad = true 
) const
virtual

Checks if a local point coordinate is inside the element.

Parameters
locpoint in local coordinates
padApply padding to the element to ensure that boundary points are considered inside.
Returns
true if the point is inside the element, false otherwise.
Note
If pad == false, then the status of points exactly on the element boundary is undefined. They may or may not be considered inside. Use pad = true to ensure that they are regarded as inside.
See Also
GContains

Implements Element.

const RVector& Triangle3old::LNormal ( int  side) const
virtual

Returns a side normal in local coordinates.

Parameters
sideside index (>= 0)
Returns
Normal to side in local element coordinates.
See Also
DirectionCosine

Implements Element.

Point Triangle3old::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

Implements Element.

RDenseMatrix Triangle3old::LocalShapeD ( const Point loc) const
virtual

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

Parameters
locpoint coordinates in the local element system
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

Implements Element.

RVector Triangle3old::LocalShapeF ( const Point loc) const
virtual

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

Parameters
locpoint 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
GlobalShapeF, LocalShapeD, GlobalShapeD

Implements Element.

int Triangle3old::nNode ( ) const
inlinevirtual

Returns the number of nodes associated with the element.

Returns
Number of nodes

Implements Element.

Point Triangle3old::NodeLocal ( int  node) const
virtual

Returns the local coordinates of an element node.

Parameters
nodenode index (>= 0)
Returns
Node coordinates in the element's local reference system.
See Also
Local, Global

Implements Element.

int Triangle3old::nSide ( ) const
inlinevirtual

Returns the number of element sides.

Returns
Number of sides
Note
For 3-D elements, this returns the number of element faces (4 for tetrahedra, 6, for voxels, etc). For 2-D elements, it returns the number of boundary segments (3 for triangles, 4 for pixels, etc.)

Implements Element.

int Triangle3old::nSideNode ( int  side) const
inlinevirtual

Returns the number of vertices associated with a side.

Parameters
sideside index (>= 0)
Returns
Number of vertices connected to the side.
See Also
nSide, nNode, nSideNode

Implements Element.

int Triangle3old::QuadRule ( int  order,
const double **  wght,
const Point **  absc 
) const
virtual

Returns the weights and abscissae of quadrature rules over the element.

Parameters
orderquadrature order (>= 1)
wghtpointer to an array of quadrature weights
abscpointer to an array of quadrature abscissae
Returns
Number of quadrature points
Note
The required order of a quadrature rule depends on the type of integral to be performed. Examples:
  • IntF: order = 1
  • IntFF, IntD: order = 2
  • intFFF, IntFD: order = 3
  • IntDD: order = 4 etc.

Reimplemented from Element.

int Triangle3old::SideNode ( int  side,
int  node 
) const
virtual

Returns relative node index for a side vertex.

Parameters
sideside index (>= 0)
nodeside node index (>= 0)
Returns
Relative node index
Note
Returns the node index for one of the vertices of a side.
side must be in the range 0 .. nSide()-1
node must be in the range 0 .. nSideNode (side)
The return value is in the range 0 .. nNode()-1
See Also
nNode, nSide, nSideNode

Implements Element.

double Triangle3old::SideSize ( int  side,
const NodeList nlist 
) const
virtual

Returns the size of an element side.

Parameters
sideside index (>= 0)
nlistmesh node list
Returns
Size of the element face (length for 2-D elements, area for 3-D elements).
See Also
Size

Reimplemented from Element.

Point Triangle3old::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 Element.

BYTE Triangle3old::Type ( ) const
inlinevirtual

Returns an element type identifier.

Returns
Element type identifier (see Element type identifiers).

Implements Element.

References ELID_TRI3OLD.


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