GemaCoreLib
The GeMA Core library
Public Member Functions | Protected Member Functions | List of all members
GmTri3Shape Class Reference

GmShape specialization for a Triangle with 3 nodes object. More...

#include <gmTriShape.h>

Inheritance diagram for GmTri3Shape:
Inheritance graph
[legend]
Collaboration diagram for GmTri3Shape:
Collaboration graph
[legend]

Public Member Functions

virtual GmCellType elemType () const
 Returns the type of this element.
 
virtual int numFunctions () const
 Returns the number of shape functions of this element type (equal to the number of nodes)
 
virtual void nodeNaturalCoord (int node, GmVector &coord) const
 Fills the coord vector with the set of natural coordinates for the reference shape function node. Node should be a value between 0 and numFunctions(). This is generally equal to the number of nodes but can be different for interface elements. More...
 
virtual void shapeValues (const GmVector &ncoord, GmVector &N) const
 Returns the shape function evaluated at (Zeta1, Zeta2, Zeta3). More...
 
virtual void shapePartials (const GmVector &ncoord, GmMatrix &dN, bool transposed=false) const
 Returns shape function partial derivatives with respect to Zeta1, Zeta2 and Zeta3, evaluated at the point (Zeta1, Zeta2, Zeta3). More...
 
virtual double edgeScalingFactor (int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
 Returns the scaling factor needed when calculating integrals over element edges. This has the same effect as the jacobian determinant when integrating over the whole element. More...
 
- Public Member Functions inherited from GmTriShape
virtual int numNaturalCoord () const
 Returns the number of natural coordinates used by this element type.
 
virtual int numCartesianCoord () const
 Returns the number of cartesian coordinates expected by this element type.
 
virtual void naturalCoordLimits (int coord, double *min, double *max) const
 Fills min and max with the domain limits for the given natural coordinate (between 0 and numNaturalCoord() - 1)
 
virtual void naturalCenter (GmVector &coord) const
 Fills the coord vector with the set of natural coordinates for the element center. More...
 
virtual bool translateEdgePoint (int edge, const GmVector &srcEdgeCoord, GmVector &elementCoord) const
 Given a "bar like" edge coordinate (from -1 to 1), at the given edge, fills elementCoord with the equivalent element natural coordinate. Should be implemented for 2D & 3D elements. More...
 
virtual bool translateFacePoint (int face, const GmVector &srcFaceCoord, GmVector &elementCoord) const
 Given a "quad like" or "tri like" face coordinate (-1 to 1 pair for quad faces and 0 to 1 barycentric triple for tri faces), at the given face, fills elementCoord with the equivalent element natural coordinate. Should be implemented for 3D elements. More...
 
virtual double shapeCartesianPartialsFromPJ (const GmMatrix &P, const GmMatrix &J, GmMatrix &dN, bool transposed=false) const
 Alternative version of shapeCartesianPartialsXxxx() to calculate shape function partial derivatives with respect to its cartesian coordinates, given a previously calculated Jacobian matrix and a previously calculated shape function partial matrix. More...
 
virtual double scaledJacobianDet (const GmMatrix &J) const
 Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differential area element. dOmega = dx.dy = J.dXi.dEta where J is the scaled Jacobian determinant. More...
 
virtual void jacobianAndPartials (const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed=false) const
 This function does the same calculations as the jacobian() call, but also filling the extra parameter P with the matrix of shape function partials. More...
 
virtual double borderScalingFactor (int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
 Returns the scaling factor needed when calculating integrals over borders (edges or faces). This has the same effect as the jacobian determinant when integrating over the whole element. It can be used to treat 2D and 3D edge / face boundary conditions in a single code. More...
 
virtual double faceScalingFactor (int border, const GmVector &borderCoord, const GmVector &elementCoord, const GmMatrix &X, bool transposed=false) const
 Returns the scaling factor needed when calculating integrals over element faces. This has the same effect as the jacobian determinant when integrating over the whole element. More...
 
- Public Member Functions inherited from GmShape
void fillNaturalCoordinates (GmMatrix &coord, bool transposed=false) const
 Fills the coord matrix with the set of natural coordinates for the shape nodes. More...
 
virtual bool cartesianToNatural (const GmVector &coord, const GmMatrix &X, GmVector &ncoord, bool *inside) const
 Given a set of cartesian coordinates 'coord' and a matrix with node coordinates 'X', calculates the corresponding natural coordinate filling 'ncoord'. The inside parameter is filled with true if the given cartesian coordinate is inside the element, false otherwise. More...
 
virtual bool cartesianToNatural (const GmVector &coord, const GmMatrix &X, bool preferGeometric, double gradTol, int maxGradIter, double outsideNatTol, GmVector &ncoord, bool *inside) const
 Given a set of cartesian coordinates 'coord' and a matrix with node coordinates 'X', calculates the corresponding natural coordinate filling 'ncoord'. The inside parameter is filled with true if the given cartesian coordinate is inside the element, false otherwise. More...
 
virtual void naturalToCartesian (const GmVector &ncoord, const GmMatrix &X, GmVector &coord) const
 Given a set of natural coordinates 'ncoord' and a matrix with node coordinates 'X', calculates the corresponding cartesian coordinates filling 'coord'. More...
 
virtual double shapeCartesianPartialsFromCoord (const GmVector &ncoord, const GmMatrix &X, GmMatrix &dN, bool transposed=false) const
 Function used to calculate shape function partial derivatives with respect to its cartesian coordinates, evaluated at the given point. More...
 
virtual double shapeCartesianPartialsFromJacobian (const GmVector &ncoord, const GmMatrix &J, GmMatrix &dN, bool transposed=false) const
 Alternative version of shapeCartesianPartials() to calculate shape function partial derivatives with respect to its cartesian coordinates, given a previously calculated Jacobian matrix. More...
 
virtual void jacobian (const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, bool transposed=false) const
 Calculates the jacobian matrix, relating cartesian coordinates to natural coordinates. More...
 
virtual const GmMatrixgaussExtrapolationMatrix (const GmIntegrationRule *ir) const
 Given an integration rule, compatible with the current element type, constructs an extrapolation matrix to calculate node values from integration point values. More...
 
virtual bool pointExtrapolationMatrix (const GmMatrix &points, GmMatrix &em) const
 Given the set of 'm' natural coordinates in points, builds an 'n x m' extrapolation matrix that can be used to extrapolate values from an 'm' sized value vector to the 'n' elements nodes. More...
 
double interpolate (const GmVector &nodeValues, const GmVector &ncoord) const
 Given a set of node values and a set of natural coordinates, interpolates the node values on the given coordinate. More...
 

Protected Member Functions

virtual bool hasGeometryBasedCartesianToNatural () const
 A virtual function that should be replaced to return true if the shape implements a geometry based algorithm for converting cartesian coordinates to natural coordinates. In that case, geometryBasedCartesianToNatural() must also be reimplemented.
 
virtual bool geometryBasedCartesianToNatural (const GmVector &coord, const GmMatrix &X, double natTol, GmVector &ncoord, bool *inside) const
 Virtual function that should be implemented if a shape can provide a geometric algorithm for converting cartesian coordinates to natural coordinates. See the cartesianToNatural() description for the arguments description.
 
- Protected Member Functions inherited from GmTriShape
virtual bool gradientBasedCartesianToNatural (const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
 
- Protected Member Functions inherited from GmShape
virtual void jacIndependentNatCoord (const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed) const
 A generic implementation of jacobianAndPartials() for the case where all natural coordinates are independent, such as in line segments, quads and hexahedrons. More...
 
virtual void jacDependentNatCoord (const GmVector &ncoord, const GmMatrix &X, GmMatrix &J, GmMatrix &P, bool transposed) const
 A generic implementation of jacobianAndPartials() for the case where all natural coordinates are NOT independent, such as in triangles and tetrahedrons. More...
 
bool ginv (const GmMatrix &A, GmMatrix &inv) const
 Calculates the Moore-Penrose generalized matrix inverse for the 'm x n' matrix A, filling inv. More...
 
bool gradientBasedCartesianToNaturalQuadHex (const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
 Given a set of cartesian coordinates 'coord' and a matrix with node coordinates 'X', calculates the corresponding natural coordinate filling 'ncoord'. Returns true if the calculation succeeded, false otherwise. More...
 
bool gradientBasedCartesianToNaturalTriTet (const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
 An implementation of the gradient based cartesian to natural algorithm specific to Tri and Tet elements. For a description of the algorithm and function parameters, see gradientBasedCartesianToNaturalQuadHex().
 
bool gradientBasedCartesianToNaturalBar (const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
 An implementation of the gradient based cartesian to natural algorithm specific to Bar elements. For a description of the algorithm and function parameters, see gradientBasedCartesianToNaturalQuadHex().
 
bool gradientBasedCartesianToNaturalWedge (const GmVector &coord, const GmMatrix &X, double tol, int maxIter, double natTol, GmVector &ncoord, bool *inside) const
 An implementation of the gradient based cartesian to natural algorithm specific to Wedge elements. For a description of the algorithm and function parameters, see gradientBasedCartesianToNaturalQuadHex().
 
bool checkAndClampNaturalCoordinate (double &ncoord, double min, double max, double natTol) const
 An utilitary function that given a natural coordinate component and its domain (given by min and max), checks that coord is inside the interval, possibly clamping it to the domain if outside by a little bit, controled by the given tolerance. More...
 
bool checkAndClampBarycentricNaturalCoordinates (GmVector &ncoord, double natTol) const
 An utilitary function that given a set of 3 or 4 barycentric natural coordinates, checks for out of the triangle/tetrahedron values, possibly clamping it to the domain if outside by a little bit, controled by the given tolerance. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from GmTriShape
static int fixedEdgeCoord (int border)
 Given a border index referencing a Tri element edge, returns the index of the natural coordinate fixed over that edge. More...
 
- Static Public Member Functions inherited from GmShape
static bool initShapeFunctions ()
 Initializes the shape function module. Must be called BEFORE any calls to shapeFromElementType() or linearShapeFromElementType(), but AFTER the static initialization of the geometry objects. More...
 
static void setConfigOptions (const GmSimulationData *simData)
 Updates the default config options used by the "simple" cartesianToNatural() call.
 
static const GmShapeshapeFromElementType (GmCellType etype, int P, int Q)
 Returns a shape object suitable for handling elements of type etype. Parameters P and Q are needed only when working with hierarchical element types, such as GM_HQUADP or GM_HHEXP. They can be set to 0 for "common" element types.
 
static const GmShapelinearShapeFromElementType (GmCellType etype, int P, int Q)
 Returns a linear shape object suitable for handling elements of type etype. Parameters P and Q are needed only when working with hierarchical element types, such as GM_HQUADP or GM_HHEXP. They can be set to 0 for "common" element types.
 

Detailed Description

GmShape specialization for a Triangle with 3 nodes object.

Natural coordinates: Zeta1, Zeta2, Zeta3 (Zeta1 + Zeta2 + Zeta3 = 1) Natural coordinates domain: [0, 1], [0, 1], [0, 1]

Node ordering: Node list formed by the three nodes defining the triangle in CCW order. See schema in the documentation of the GmCellGeometry class.

Member Function Documentation

◆ edgeScalingFactor()

double GmTri3Shape::edgeScalingFactor ( int  border,
const GmVector borderCoord,
const GmVector elementCoord,
const GmMatrix X,
bool  transposed = false 
) const
virtual

Returns the scaling factor needed when calculating integrals over element edges. This has the same effect as the jacobian determinant when integrating over the whole element.

This function is undefined for line elements.

Parameters
borderDefines the edge where the integration is taking place.
borderCoordDefines the integration point where the factor will be calculated in the edge reference system.
elementCoordDefines the integration point where the factor will be calculated in the element reference system. It should have a dimension equal to numNaturalCoord().
XThe X matrix gives the cartesian coordinates of every node and should be organized as a 'd x n' matrix wher d is the number of dimension and n the number of shape functions. See the description given in jacobian().
transposedIf true, the X matrix should be the transposed of the matrix defined above.
Returns
Returns the calculated scaling factor

Implements GmShape.

◆ nodeNaturalCoord()

void GmTri3Shape::nodeNaturalCoord ( int  node,
GmVector coord 
) const
virtual

Fills the coord vector with the set of natural coordinates for the reference shape function node. Node should be a value between 0 and numFunctions(). This is generally equal to the number of nodes but can be different for interface elements.

The resulting vector will have size equal to numNaturalCoord(). If coord holds a vector created by DECLARE_REF_VECTOR(), its size MUST already be equal to the expected result size.

Implements GmShape.

◆ shapePartials()

void GmTri3Shape::shapePartials ( const GmVector ncoord,
GmMatrix dN,
bool  transposed = false 
) const
virtual

Returns shape function partial derivatives with respect to Zeta1, Zeta2 and Zeta3, evaluated at the point (Zeta1, Zeta2, Zeta3).

See comments on the base class for function behaviour and on the class for coordinate assumptions and node ordering

Tri 3 shape functions (Felippa, IFEM equation 16.9, page 16-6):

N1 = Z1 N2 = Z2 N3 = Z3

Partial derivatives with respect to Zeta1: dN1/dZ1 = 1 dN2/dZ1 = 0 dN3/dZ1 = 0

Partial derivatives with respect to Zeta2: dN1/dZ2 = 0 dN2/dZ2 = 1 dN3/dZ2 = 0

Partial derivatives with respect to Zeta3: dN1/dZ3 = 0 dN2/dZ3 = 0 dN3/dZ3 = 1

Implements GmShape.

◆ shapeValues()

void GmTri3Shape::shapeValues ( const GmVector ncoord,
GmVector N 
) const
virtual

Returns the shape function evaluated at (Zeta1, Zeta2, Zeta3).

See comments on the base class for function behaviour and on the class for coordinate assumptions and node ordering

Tri 3 shape functions (Felippa, IFEM equation 16.9, page 16-6):

N1 = Zeta1 N2 = Zeta2 N3 = Zeta3

Implements GmShape.


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