GeMA
The GeMA main application
Shape function object methods

Shape function objects are a wrapper over the GeMA library providing shape function definitions for all the available Element types. Besides calculating shape function values, it also provides methods for returning their partial derivatives, with respect to both natural and cartesian coordinates, as well as methods for calculating Jacobian matrices.

Shape function objects are associated with an element type, and not with a particular element. A shape function object for a given element type can be retrieved by calling the ShapeFunction() function. Alternatively, the shape function object associated with an element can be retrieved by calling the element:shape() and element:linearShape() methods.

Example:

-- Get the shape function object for the 'hex8' element type
local shape = ShapeFunction('hex8')
-- Get a reference to the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e'
local shape = e:shape()

Most used:

Index:

Creating a shape function object

ShapeFunction(elemType, linear)
ShapeFunction(elemType, P, Q, linear)
Description: Returns the shape function object associated with the given element type.
Parameters: elemType A string with the element type ("quad4", "tri3", ...). See the Element types page for a list of available element types.
P, Q The requested element polynomial interpolation order when elemType is an hierarchical element. P gives the order used for scalar fields and Q for vector fields.
linear An optional flag. When elemType refers to a quadratic element, and this flag is set to true, returns the shape function from the equivalent linear element (for a quad9 element, for example, returns a quad4 shape function). Usefull when working with super parametric element types. Default = false.
Returns: Returns the shape function object.

Example:

local shape = ShapeFunction('quad8') -- quad8 shape function
local lshape = ShapeFunction('quad8', true) -- quad4 shape function
local hshape = ShapeFunction('hquadp', 3, 2) -- hquadP shape function with interpolation orders P = 3 and Q = 2


Shape function metadata methods

shape:elemType()
Description: Returns the associated element type.
Parameters: None.
Returns: Returns a string with the associated element type ("quad4", "tri3", ...). See the Element types page for a list of available types.

Example:

local type = shape:elemType()


shape:numFunctions()
Description: Returns the number os shape functions associated with this element type. This value is usually equal to the number of associated element nodes, but can be different for interface elements.
Parameters: None.
Returns: Returns the number of shape functions for this type.

Example:

local nfunctions = shape:numFunctions()


shape:numNaturalCoord()
Description: Returns the dimension for this shape function natural coordinates. This value might be different from the number of Cartesian coordinates. Shape functions for a tri3 element, for example, use 3 natural coordinates against 2 Cartesian coordinates.
Parameters: None.
Returns: Returns the number of natural coordinates used by this element type.

Example:

local nnatural = shape:numNaturalCoord()


shape:numCartesianCoord()
Description: Returns the expected dimension for node Cartesian coordinates. Shape functions for an hex8 element, for example, expect node coordinates with 3 dimensions.
Parameters: None.
Returns: Returns the number of Cartesian coordinates expected by this element type.

Example:

local ncartesian = shape:numCartesianCoord()


shape:naturalCoordLimits(coordDim)
Description: Returns the domain limits (min, max) for the given natural coordinate dimension. Each element type can work in its preferred domain. A quad4 element, for example, works with a [-1, 1] domain, while a tri3 element works with a [0, 1] domain. The coordinate dimension parameter is necessary since elements like a wedge might use different domains for each of its natural coordinate components.
Parameters: coordDim - Coordinate dimension index (between 1 and shape:numNaturalCoord()).
Returns: Returns two numbers with the domain interval for the given coordinate dimension.

Example:

local min, max = shape:naturalCoordLimits(1)


Coordinate retrieval and conversion methods

shape:nodeNaturalCoord(node)
Description: Returns a vector object filled with the natural coordinates for the given local element node index.
Parameters: node - A local node index. In reallity this parameter is really a shape function index (between 1 and shape:numFunctions()) since in some cases (interface elements) the number of element nodes is not equal to the number of shape functions.
Returns: Returns a vector object with the node's natural coordinates. Its size will be equal to shape:numNaturalCoord().

Example:

-- Get the natural coordinates for the first element node
local ncoord = shape:nodeNaturalCoord(1)


shape:naturalCenter()
Description: Returns a vector object filled with the natural coordinates for the element center.
Parameters: None.
Returns: Returns a vector object with the element's center natural coordinates. Its size will be equal to shape:numNaturalCoord().

Example:

-- Get the natural coordinates for the element center
local ncoord = shape:naturalCenter()


shape:naturalCoordinates(transposed)
Description: Returns a matrix object filled with the natural coordinates of the element's nodes.
Parameters: transposed - Optional parameter specifying that the result matrix should be transposed.
Returns: Returns a n x m matrix object with node natural coordinates, where n is the number of shape functions (usually equal to the number of element nodes - see comments on shape:nodeNaturalCoord()) and m the number of natural coordinates.
Returns a m x n matrix if transposed is true.
E.g.: quad4, n=4, m=2: \(\small\begin{pmatrix} -1.0 & -1.0 \\ 1.0 & -1.0 \\ 1.0 & 1.0 \\ -1.0 & 1.0 \end{pmatrix}\)

Example:

-- Get the natural coordinates for every element node
local ncoordMatrix = shape:naturalCoordinates()


shape:translateEdgePoint(edge, srcEdgeCoord)
Description: Given an edge and a "bar like" edge coordinate, a value between -1 and 1, returns that coordinate translated to the element natural coordinate reference system.
Parameters: edge The edge index, a value between 1 and the number of element edges (geometry:numEdges()). Check the Element types page for edge numbering for each cell type.
srcEdgeCoord A vector object with the edge natural coordinate in a "bar like" system (a value from -1 to 1).
Returns: Returns the vector object with the given point translated to the element natural coordinate system or nil if edge references an element edge that can not be used for numerical integration (like the 0-length edges for an interface element).

Example:

local quadCoord = ShapeFunction('quad4'):translateEdgePoint(2, Vector{0.5})
assert(quadCoord == Vector{1, 0.5})


shape:translateFacePoint(face, srcFaceCoord)
Description: Given a face and a "quad or tri like" face coordinate, a value pair between -1 and 1 for quad faces and a barycentric triple coordinate, between 0 and 1, for triangular faces, returns that coordinate translated to the element natural coordinate reference system.
Parameters: face The face index, a value between 1 and the number of element faces (geometry:numFaces()). Check the Element types page for face numbering for each cell type.
srcFaceCoord A vector object with the face natural coordinate in a "quad or tri like" system (a pair of values from -1 to 1 for quadrilateral faces or a triple of barycentric coordinates for triangular faces).
Returns: Returns the vector object with the given point translated to the element natural coordinate system or nil if face references an element face that can not be used for numerical integration (like the 0-area faces for an interface element).

Example:

local hexCoord = ShapeFunction('hex8'):translateFacePoint(5, Vector{-0.5, 0.5})
assert(hexCoord == Vector{-1, -0.5, 0.5})


shape:cartesianToNatural(coord, coordMatrix)
shape:cartesianToNatural(coord, coordMatrix, preferGeometric, gradTol, maxGradIter, outTol)
Description: Given a point described by its Cartesian coordinates and a matrix with node coordinates for an element, returns a boolean stating if the point is inside the element and a vector object with the point's natural coordinates. When using the first, simplified, signature, the missing arguments from the extended signature are filled with default values specified either in the configuration file or in the model options (see the defCtoNPreferGeometric, defCtoNGradTol, defCtoNMaxGradIter and defCtoNOutsideNatTol option descriptions). Important: when using the extended signature, all the extra parameters should be supplied in the call.
Parameters: coord A vector object with the Cartesian coordinates of the point.
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of nodes in the element. Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
preferGeometric When using the extended call signature, defines the algorithm that will be used by the call. When set to true, if the element type supports a geometric based algorithm, that method will be used. Otherwise, a generic gradient based algorithm will be used. When set to false, will always use the gradient based algorithm. When using the simplified call signature, the default value given by the defCtoNPreferGeometric option is assumed.
gradTol When using the extended call signature, defines the convergence tolerance used by the gradient based algorithm. Unused for geometric algorithms. When using the simplified call signature, the default value given by the defCtoNGradTol option is assumed.
maxGradIter When using the extended call signature, defines the maximum number of iterations used by the gradient based algorithm. Unused for geometric algorithms. When using the simplified call signature, the default value given by the defCtoNMaxGradIter option is assumed.
outTol When using the extended call signature, defines a tolerance that will be used when treating coordinates outside, but very close to the element border. When greater than 0.0, this value gives a "snap" tolerance within which cartesian coordinates slightly outside the element will be snapped to the element border when converting cartesian to natural coordinates. For example, if the tolerance is 1e-3 and the conversion yields a natural coordinate of 1.0002, the coordinate will be changed to 1.0 and the point considered internal to the element. When set to 0.0, no snapping will be done. When using the simplified call signature, the default value given by the defCtoNOutsideNatTol option is assumed.
Returns: Returns true or false depending if the point is inside the element or not. Also returns a vector object with the point's natural coordinate inside the element. If the natural coordinate could not be calculated (no convergance, for example), an error is raised, aborting the simulation.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e'
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Check if the point {12.75, 7.5} is inside the first mesh element and get
-- its natural coordinates inside that element
local inside, ncoord = shape:cartesianToNatural(Vector{12.5, 7.5}, X)
-- Check if the point {12.75, 7.5} is inside the first mesh elemennt, forcing the use
-- of the gradient based algorithm and with a "snap" tolerance of 0.001
local inside, ncoord = shape:cartesianToNatural(Vector{12.5, 7.5}, X, false, 1e-8, 20, 0.001)


shape:naturalToCartesian(ncoord, coordMatrix)
Description: Given a point described by its natural coordinates and a matrix with node coordinates for the element that contains the given point, returns the point's corresponding Cartesian coordinates.
Parameters: ncoord A vector object with the natural coordinates of the point.
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of nodes in the element. Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
Returns: Returns a vector object with the point's Cartesian coordinates.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e'
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the Cartesian coordinates for the centroid of the first mesh element
local coord = shape:naturalToCartesian(Vector{0.0, 0.0}, X)


Shape values, partials and Jacobian

shape:shapeValues(ncoord)
Description: Returns the set of shape functions evaluated over a point defined by its natural coordinates.
Parameters: ncoord - A vector object with the natural coordinates of the evaluation point.
Returns: Returns a vector object filled with the calculated shape function values. Its size is equal to shape:numFunctions().

Example:

local quad4shape = ShapeFunction('quad4')
local tri3shape = ShapeFunction('tri3')
-- Evaluates the shape functions at the element centroid
local Nquad = quad4shape:shapeValues(Vector{0.0, 0.0})
local Ntri = tri3shape:shapeValues(Vector{1/3, 1/3, 1/3})


shape:shapePartials(ncoord, transposed)
Description: Returns the set of shape function partial derivatives with respect to its natural coordinates, evaluated over a point defined by its natural coordinates.
Parameters: ncoord A vector object with the natural coordinates of the evaluation point.
transposed Optional parameter specifying that the result matrix should be transposed.
Returns: Returns a n x m matrix object filled with the calculated shape function partial derivatives, where n is the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()) and m is the number of natural coordinates. In this way, line i holds the set of partial derivatives of the shape function Ni with respect to all of its natural coordinates.
Returns a m x n matrix if transposed is true.
E.g.: quad4, n=4, m=2: \(\small\begin{pmatrix} dN_1/d\xi & dN_1/d\eta \\ dN_2/d\xi & dN_2/d\eta \\ dN_3/d\xi & dN_3/d\eta \\ dN_4/d\xi & dN_4/d\eta \end{pmatrix}\)

Example:

-- Evaluates the partial derivatives of a hexahedron at its centroid
local shape = ShapeFunction('hex8')
local dN = shape:shapePartials(Vector{0.0, 0.0, 0.0})


shape:jacobian(ncoord, coordMatrix, transposed)
Description: Retuns the shape function Jacobian matrix, relating Cartesian coordinates to natural coordinates, evaluated over a point defined by its natural coordinates.
Parameters: ncoord A vector object with the natural coordinates of the evaluation point.
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()). Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
transposed Optional parameter specifying that the result matrix should be transposed. If set to true, the matrix in coordMatrix should also be transposed (becoming a n x d matrix).
Returns: Returns a n x n matrix object filled with the calculated Jacobian matrix, where n is the number of natural coordinates, organized with gradients in lines. For element types where the number of cartesian coordinates is different from the number of natural coordinates (triangle and tetrahedron families of elements), the first matrix line will be filled with ones. See Felippa, Advanced FEM, Appendix I, table I.3 for more details.
If transposed is true, the returned matrix will be transposed (gradients in columns).
E.g.: quad4, n=2: \(\begin{pmatrix} dx/d\xi & dx/d\eta \\ dy/d\xi & dy/d\eta \end{pmatrix}\)

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (a quad4)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the Jacobian matrix evaluated at the centroid of the first mesh element
local J = shape:jacobian(Vector{0.0, 0.0}, X)


shape:jacobianAndPartials(ncoord, coordMatrix, transposed)
Description: Returns the shape function Jacobian matrix and the shape function partial derivatives with respect to natural coordinates, evaluated over a point defined by its natural coordinates. Equivalent to calling shape:jacobian() and shape:shapePartials(). As both are commonly used together, this function saves the need to compute the partials matrix twice (since it is used in the Jacobian computation).
Parameters: See the description for shape:jacobian().
Returns: Returns the Jacobian matrix followed by the partials matrix. See the documentation for shape:jacobian() and shape:shapePartials() for a description of both matrix formats.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (a quad4)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the Jacobian matrix and the shape partials matrix evaluated at the centroid of the first mesh element
local J, P = shape:jacobianAndPartials(Vector{0.0, 0.0}, X)


shape:scaledJacobianDet(J)
Description: Returns the Jacobian determinant multiplied by the scaling factor needed for transforming the differential area term. \(d\Omega = dx dy = J d\xi d\eta\), where J is the value returned by this function. See Felippa Advanced FEM, appendix I, table I.3 / section I.2.4 or Felippa Introduction to FEM, chapters 17 (eq 17.20) and 24 (eq 24.7).
Parameters: J - The jacobian matrix object returned by shape:jacobian() or by shape:jacobianAndPartials().
Returns: Returns the Jacobian matrix determinant, scaled by the apropriate constant.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (a quad4)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the Jacobian matrix evaluated at the centroid of the first mesh element
-- and its scaled determinant
local J = shape:jacobian(Vector{0.0, 0.0}, X)
local detJ = shape:scaledJacobianDet(J)


shape:borderScalingFactor(border, borderCoord, elementCoord, coordMatrix, transposed)
Description: Returns the scaling factor needed for transforming the differential area term when calculating integrals over borders (edges or faces) at the supplied natural coordinate. The returned value has the same purpose as the value returned by shape:scaledJacobianDet() when integrating over the entire element. For 2D elements it is equivalent to calling shape:edgeScalingFactor() and for 3D elements to calling shape:faceScalingFactor(), beeing undefined for line elements. It can be used for building generic code that works for both 2D and 3D elements.
Parameters: border The border index, a value between 1 and the number of element edges (geometry:numEdges()) for 2D elements, and between 1 and the number of element faces (geometry:numFaces()) for 3D elements. Check the Element types page for edge and face numberings for each cell type.
borderCoord A vector object with the integration point where the factor will be calculated in the border (edge or face) natural coordinates system.
elementCoord A vector object with the integration point where the factor will be calculated in the element natural coordinate system. It should have a dimension equal to shape:numNaturalCoord(). Also, the value at the position returned by shape:fixedBorderCoord(border) should be equal to shape:fixedBorderValue(border).
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()). Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
transposed Optional parameter specifying if the matrix in coordMatrix should be a transposed version of the default version described above, thus becoming a n x d matrix.
Returns: Returns the calculated scaling factor.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (a quad4)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the scaling factor for a point evaluated at the middle of the first edge of the first mesh element.
local factor = shape:borderScalingFactor(1, Vector{0.0}, Vector{0.0, -1.0}, X)


shape:faceScalingFactor(border, borderCoord, elementCoord, coordMatrix, transposed)
Description: Returns the scaling factor needed for transforming the differential area term when calculating integrals over faces at the supplied natural coordinate. The returned value has the same purpose as the value returned by shape:scaledJacobianDet() when integrating over the entire element. Undefined for line and surface elements.
Parameters: border The border index, a value between 1 and the number of element faces (geometry:numFaces()). Check the Element types page for face numbering for each cell type.
borderCoord A vector object with the integration point where the factor will be calculated in the face natural coordinates system.
elementCoord A vector object with the integration point where the factor will be calculated in the element natural coordinate system. It should have a dimension equal to shape:numNaturalCoord(). Also, the value at the position returned by shape:fixedFaceCoord(border) should be equal to shape:fixedFaceValue(border).
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()). Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
transposed Optional parameter specifying if the matrix in coordMatrix should be a transposed version of the default version described above, thus becoming a n x d matrix.
Returns: Returns the calculated scaling factor.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (an hex8)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the scaling factor for a point evaluated at the middle of the first face of the first mesh element.
local factor = shape:faceScalingFactor(1, Vector{0.0, 0.0}, Vector{0.0, 0.0, -1.0}, X)


shape:edgeScalingFactor(border, borderCoord, elementCoord, coordMatrix, transposed)
Description: Returns the scaling factor needed for transforming the differential area term when calculating integrals over edges at the supplied natural coordinate. The returned value has the same purpose as the value returned by shape:scaledJacobianDet() when integrating over the entire element. Undefined for line elements.
Parameters: border The edge index, a value between 1 and the number of element edges (geometry:numEdges()). Check the Element types page for edge numbering for each cell type.
borderCoord A vector object with the integration point where the factor will be calculated in the edge natural coordinates system.
elementCoord A vector object with the integration point where the factor will be calculated in the element natural coordinate system. It should have a dimension equal to shape:numNaturalCoord(). Also, the value(s) at the position(s) returned by shape:fixedEdgeCoord(border) should be equal to shape:fixedEdgeValue(border).
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()). Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
transposed Optional parameter specifying if the matrix in coordMatrix should be a transposed version of the default version described above, thus becoming a n x d matrix.
Returns: Returns the calculated scaling factor.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (a quad4)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the scaling factor for a point evaluated at the middle of the first edge of the first mesh element.
local factor = shape:edgeScalingFactor(1, Vector{0.0}, Vector{0.0, -1.0}, X)


shape:shapeCartesianPartialsFromCoord(ncoord, coordMatrix, transposed)
Description: Returns the set of shape function partial derivatives with respect to its Cartesian coordinates, evaluated over a point defined by its natural coordinates. Also returns the scaled Jacobian determinant.
Parameters: ncoord A vector object with the natural coordinates of the evaluation point.
coordMatrix A matrix object with node coordinates in a d x n format, where d is the dimension of the Cartesian coordinates and n the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()). Can be obtained by a call to elem:nodeMatrix(accessor, true). Notice that a true value must be passed as second parameter so that the returned matrix conforms to the expected layout. E.g.: quad4, d=2, n=4: \(\begin{pmatrix} x_1 & x_2 & x_3 & x_4 \\ y_1 & y_2 & y_3 & y_4 \end{pmatrix}\)
transposed Optional parameter specifying that the result matrix should be transposed. If set to true, the matrix in coordMatrix should also be transposed (becoming a n x d matrix).
Returns: Returns a n x d matrix object filled with the calculated shape function Cartesian partial derivatives, where n is the number of shape functions (usually equal to the number of element nodes - see comments on shape:numFunctions()) and d is the number of Cartesian coordinates. In this way, line i holds the set of partial derivatives of the shape function Ni with respect to all of its Cartesian coordinates.
Returns a d x n matrix if transposed is true.
This function also returns the scaled Jacobian determinant, often used together with the shape Cartesian partials. This is the same value returned by a call to shape:scaledJacobianDet() but reusing the Jacobian matrix used while calculating the shape Cartesian partials.
E.g.: quad4, n=4, d=2: \(\small\begin{pmatrix} dN_1/dx & dN_1/dy \\ dN_2/dx & dN_2/dy \\ dN_3/dx & dN_3/dy \\ dN_4/dx & dN_4/dy \end{pmatrix}\)

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the shape function object associated with element 'e' (a quad4)
local shape = e:shape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the shape Cartesian partials matrix and the Jacobian determinant evaluated at the centroid of the first mesh element
local B, detJ = shape:shapeCartesianPartialsFromCoord(Vector{0.0, 0.0}, X)


shape:shapeCartesianPartialsFromJacobian(ncoord, J, transposed)
Description: Returns the set of shape function partial derivatives with respect to its Cartesian coordinates, evaluated over a point defined by its natural coordinates, given a previously calculated Jacobian matrix. Also returns the scaled Jacobian determinant.
This function is usefull when working with super-parametric elements where the Jacobian is calculated using the 'native' element shape functions but field values are calculated with a linear shape function.
Parameters: ncoord A vector object with the natural coordinates of the evaluation point.
J The Jacobian matrix, calculated by a call to shape:jacobian(), passing as argument the same transposed parameter received by this function.
transposed Optional parameter specifying that the result matrix should be transposed. If set to true, the Jacobian matrix in J should also be transposed.
Returns: Returns the Cartesian partial derivatives matrix and the scaled Jacobian determinant. See the shape:shapeCartesianPartialsFromCoord() method for further details.

Example:

-- Get a reference for the first mesh element
local m = modelData:mesh('myMeshName')
local e = m:cell(1)
-- Get the element and linear shape functions associated with element 'e' (a quad8)
local elemShape = e:shape()
local linearShape = e:linearShape()
-- Get the node coordinate accessor
local coordAc = m:nodeCoordAccessor()
-- Get the matrix with node coordinates
local X = e:nodeMatrix(coordAc, true)
-- Get the shape Cartesian partials matrix and the Jacobian determinant evaluated at the centroid
-- of the first mesh element. Since this is a super-parametric scenario, calculate sthe Jacobian with
-- the complete set of shape functions, and the B matrix with the linear one.
local J = elemShape:jacobian(Vector{0.0, 0.0}, X)
local B, detJ = linearShape:shapeCartesianPartialsFromJacobian(Vector{0.0, 0.0}, J)


Interpolation and extrapolation methods

shape:interpolate(values, ncoord)
Description: Given a set of node scalar values and the natural coordinates of the interpolation point, calculates the value on the given coordinate by using the shape functions to interpolate the given nodal values. For additional ways for interpolating values inside an element from nodal data, including multi-dimensional data, please refer to the Interpolation functions page.
Parameters: values A vector object with the given node values. Should have size equal to the number of element nodes (cell:numNodes()) and must be organized following the element node order convention (see the Element types page).
Parameters: ncoord A vector object with the interpolation point natural coordinates (size equal to shape:numNaturalCoord()).
Returns: Returns the scalar interpolated value.

Example:

-- Interpolates the value on the centroid of a quad element given that nodal values are 1, 5, 7 and 3.
local shape = ShapeFunction('quad4')
local v = shape:interpolate(Vector{1, 5, 7, 3}, Vector{0, 0})


shape:gaussExtrapolationMatrix(ir)
Description: Returns an extrapolation matrix appropriate for extrapolating values on an element's Gauss points to the element's nodes. For that, the resulting matrix must be multiplied by a vector filled with the values at the Gauss points, ordered following the iteration order defined by the ir:integrationPoint() method. For additional ways for extrapolating Gauss values to nodes, please refer to the Interpolation functions page.
Parameters: ir - The element integration rule object specifying the set of Gauss points where values are known.
Returns: Returns the calculated n x m extrapolation matrix, where n is the number of nodes in the element and m is the number of Gauss points in the given rule, or an empty matrix if this function is not supported for this element type / integration rule pair.

Example:

local shape = ShapeFunction('quad4')
local m = shape:gaussExtrapolationMatrix(shape:integrationRule('auto', 3, 3)) -- 3 x 3 rule
assert(m:nlin() == 4)
assert(m:ncol() == 9)
local m2 = s1:gaussExtrapolationMatrix(s1:integrationRule('auto', 1, 4)) -- Unsupported rule
assert(m2:nlin() == 0 and m2:ncol() == 0)