fem#

The EasyFEA/fem/ module in EasyFEA provides essential tools for creating and managing finite element meshes, which are crucial for numerical simulations using the Finite Element Method (FEM).


What is a mesh in EasyFEA?#

A Mesh object in EasyFEA represents a collection of ElemType used to define the geometry and structure for finite element analysis. It contains multiple _GroupElem instances, which are groups of ElemType that collectively define the spatial discretization of the domain for numerical simulations.

For example, a HEXA8 mesh includes the following element types:

All implemented element types, along with their corresponding shape functions and derivatives, are defined in the EasyFEA/fem/elems/ directory. The Gauss point quadratures are implemented in the EasyFEA/fem/_gauss.py module.


Creating or importing a Mesh#

To construct a Mesh using the Mesher, you must first create _Geom objects (see geoms for some examples). The Mesher class serves as an interface to Gmsh, a powerful meshing tool, and includes the following primary functions for mesh generation:

Several examples are available in Meshes.


Detailed fem API#

class EasyFEA.fem.BiLinearForm(form)[source]#

Bases: _Form

Assemble(field)[source]#

Assemble de form with the field.

Parameters:

field (Field) – field

Returns:

the assembled (Nn * dof_n, Nn * dof_n) sparse matrix.

Return type:

csr_matrix

Integrate_e(field)[source]#

Integrates de form with the field on elements.

Parameters:

field (Field) – field

Returns:

the integrated (Ne, …) numpy array

Return type:

np.ndarray

class EasyFEA.fem.BoundaryCondition(problemType, nodes, dofs, unknowns, dofsValues, description)[source]#

Bases: object

static Get_dofs(problemType, list_Bc_Condition)[source]#

Returns the degrees of freedom for the problem type.

Parameters:
  • problemType (str) – Problem type.

  • list_Bc_Condition (list[BoundaryCondition]) – List of boundary conditions.

Returns:

degrees of freedom.

Return type:

_types.IntArray

static Get_dofs_nodes(availableUnknowns, nodes, unknowns)[source]#

Retrieves degrees of freedom (dofs) associated with the nodes.

Parameters:
  • availableUnknowns (list[str]) – Available dofs as a list of strings. Must be a unique string list.

  • nodes (_types.IntArray) – Nodes for which dofs are calculated.

  • unknowns (list[str]) – unknowns.

Returns:

degrees of freedom.

Return type:

_types.IntArray

static Get_nBc(problemType, list_Bc_Condition)[source]#

Returns the number of conditions for the problem type.

Parameters:
  • problemType (str) – Problem type.

  • list_Bc_Condition (list[BoundaryCondition]) – List of boundary conditions.

Returns:

Number of boundary conditions (nBc).

Return type:

int

static Get_values(problemType, list_Bc_Condition)[source]#

Returns the dofs values for problem type.

Parameters:
  • problemType (str) – Problem type.

  • list_Bc_Condition (list[BoundaryCondition]) – List of boundary condition.

Returns:

dofs values.

Return type:

_types.FloatArray

property dofs: ndarray[tuple[Any, ...], dtype[IntType]]#

degrees of freedom associated with the nodes and unknowns

property dofsValues: ndarray[tuple[Any, ...], dtype[floating]]#

values applied

property nodes: ndarray[tuple[Any, ...], dtype[IntType]]#

nodes on which the condition is applied

property problemType: str#

type of problem

property unknowns: list[str]#

dofs unknowns

EasyFEA.fem.Calc_projector(oldMesh, newMesh)[source]#

Get the matrix used to project the solution from the old mesh to the new mesh such that:

newU = proj * oldU

(newNn) = (newNn x oldNn) (oldNn)

Parameters:
  • oldMesh (Mesh) – old mesh

  • newMesh (Mesh) – new mesh

Returns:

projection matrix (newMesh.Nn, oldMesh.Nn)

Return type:

sp.csr_matrix

EasyFEA.fem.Det(mat)[source]#

Computes det(mat)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

class EasyFEA.fem.ElemType(value)[source]#

Bases: str, Enum

Implemented Lagrange isoparametric element types.

static Get_1D()[source]#

Returns 1D element types.

Return type:

list[ElemType]

static Get_2D()[source]#

Returns 2D element types.

Return type:

list[ElemType]

static Get_3D()[source]#

Returns 3D element types.

Return type:

list[ElemType]

HEXA20 = 'HEXA20'#
HEXA27 = 'HEXA27'#
HEXA8 = 'HEXA8'#
POINT = 'POINT'#
PRISM15 = 'PRISM15'#
PRISM18 = 'PRISM18'#
PRISM6 = 'PRISM6'#
QUAD4 = 'QUAD4'#
QUAD8 = 'QUAD8'#
QUAD9 = 'QUAD9'#
SEG2 = 'SEG2'#
SEG3 = 'SEG3'#
SEG4 = 'SEG4'#
SEG5 = 'SEG5'#
TETRA10 = 'TETRA10'#
TETRA4 = 'TETRA4'#
TRI10 = 'TRI10'#
TRI15 = 'TRI15'#
TRI3 = 'TRI3'#
TRI6 = 'TRI6'#
property topology: str#
class EasyFEA.fem.FeArray(input_array, broadcastFeArrays=False)[source]#

Bases: ndarray[tuple[Any, …], dtype[Any]]

Finite Element array.

FeArray is a Python class designed to optimize finite element simulations by leveraging NumPy arrays with a shape of (Ne, nPg, …). This structure enables vectorized operations, eliminating the need for slow loops over elements and integration points. By using np.einsum, it efficiently handles tensor computations, significantly improving performance and code clarity for finite element analyses.

FeArrayALike#

alias of FeArray | ndarray[tuple[Any, …], dtype[Any]]

property T: FeArray | ndarray[tuple[Any, ...], dtype[Any]]#

View of the transposed array.

Same as self.transpose().

Examples

>>> import numpy as np
>>> a = np.array([[1, 2], [3, 4]])
>>> a
array([[1, 2],
       [3, 4]])
>>> a.T
array([[1, 3],
       [2, 4]])
>>> a = np.array([1, 2, 3, 4])
>>> a
array([1, 2, 3, 4])
>>> a.T
array([1, 2, 3, 4])

See also

transpose

static asfearray(array, broadcastFeArrays=False)[source]#
Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

ddot(other)[source]#
Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

dot(other)[source]#

Refer to numpy.dot() for full documentation.

See also

numpy.dot

equivalent function

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

max(*args, **kwargs)[source]#

np.max() wrapper.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

min(*args, **kwargs)[source]#

np.min() wrapper.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

static ones(*shape, dtype=None)[source]#
Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

sum(*args, **kwargs)[source]#

np.sum() wrapper.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

static zeros(*shape, dtype=None)[source]#
Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

class EasyFEA.fem.Field(groupElem, dof_n, matrixType=MatrixType.mass)[source]#

Bases: object

Field class.

Evaluate_e(function, dofsValues, returnMeanValues=True)[source]#

Evaluates the given function for the provided field over the elements.

Parameters:
  • function (Callable[[Field], FeArray]) – A function that takes a Field as input and returns a FeArray.

  • dofsValues (np.ndarray) – Array of shape (Nn * field.dof_n,) containing the degrees of freedom values.

  • returnMeanValues (bool, optional) – If True, returns the mean of the values at each element. Default is True.

Return type:

ndarray

Evaluate_n(function, dofsValues)[source]#

Evaluates the given function for the provided field over the nodes.

Parameters:
  • function (Callable[[Field], FeArray]) – A function that takes a Field as input and returns a FeArray.

  • dofsValues (np.ndarray) – Array of shape (Nn * field.dof_n,) containing the degrees of freedom values.

Return type:

ndarray

Get_coords(concatenate=False)[source]#

Returns integration point coordinates (x,y,z) for each element.

Interpolate(dofsValues)[source]#

Interpolates degrees of freedom values at each integration point for every element.

Parameters:

dofsValues (np.ndarray) – Array of shape (Nn * dof_n,) containing the degrees of freedom values.

Returns:

The (Ne, nPg, dof_n) finite element array.

Return type:

FeArray

copy()[source]#
Return type:

Field

ddot(other)[source]#
property dof_n: int#

degrees of freedom per node.

dot(other)[source]#
property grad: FeArray | ndarray[tuple[Any, ...], dtype[Any]]#

Returns the gradient of the field.

property groupElem: _GroupElem#

Group of elements.

property matrixType: MatrixType#
class EasyFEA.fem.Gauss(elemType, matrixType)[source]#

Bases: object

static Gauss_factory(elemType, matrixType)[source]#

Returns the integration points and weights based on the element and matrix type.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]

property coord: ndarray[tuple[Any, ...], dtype[floating]]#

integration point coordinates

property nPg: int#

number of integration points

property weights: ndarray[tuple[Any, ...], dtype[floating]]#

integration point weights

class EasyFEA.fem.GroupElemFactory[source]#

Bases: object

static Create(elemType, connect, coordGlob)[source]#

Creates an element group

Parameters:
  • elemType (ElemType) – element type

  • connect (_types.IntArray) – connection matrix storing nodes for each element (Ne, nPe)

  • coordGlob (_types.FloatArray) – nodes coordinates

Returns:

the element group

Return type:

GroupeElem

DICT_ELEMTYPE: dict[ElemType, tuple[int, int, int, int, int, int, int, int]] = {ElemType.HEXA20: (17, 20, 3, 2, 8, 12, 0, 0), ElemType.HEXA27: (12, 27, 3, 2, 8, 12, 6, 1), ElemType.HEXA8: (5, 8, 3, 1, 8, 0, 0, 0), ElemType.POINT: (15, 1, 0, 0, 0, 0, 0, 0), ElemType.PRISM15: (18, 15, 3, 2, 6, 9, 0, 0), ElemType.PRISM18: (13, 18, 3, 2, 6, 9, 3, 0), ElemType.PRISM6: (6, 6, 3, 1, 6, 0, 0, 0), ElemType.QUAD4: (3, 4, 2, 1, 4, 0, 0, 0), ElemType.QUAD8: (16, 8, 2, 2, 4, 4, 0, 0), ElemType.QUAD9: (10, 9, 2, 2, 4, 4, 1, 0), ElemType.SEG2: (1, 2, 1, 1, 2, 0, 0, 0), ElemType.SEG3: (8, 3, 1, 2, 2, 1, 0, 0), ElemType.SEG4: (26, 4, 1, 3, 2, 2, 0, 0), ElemType.SEG5: (27, 5, 1, 4, 2, 3, 0, 0), ElemType.TETRA10: (11, 10, 3, 2, 4, 6, 0, 0), ElemType.TETRA4: (4, 4, 3, 1, 4, 0, 0, 0), ElemType.TRI10: (21, 10, 2, 3, 3, 6, 1, 0), ElemType.TRI15: (23, 15, 2, 4, 3, 9, 3, 0), ElemType.TRI3: (2, 3, 2, 1, 3, 0, 0, 0), ElemType.TRI6: (9, 6, 2, 2, 3, 3, 0, 0)}#

(gmshId, nPe, dim, order, Nvertex, Nedge, Nface, Nvolume)

Type:

ElemType

DICT_GMSH_DATA: dict[int, tuple[ElemType, int, int, int, int, int, int, int]] = {1: (ElemType.SEG2, 2, 1, 1, 2, 0, 0, 0), 2: (ElemType.TRI3, 3, 2, 1, 3, 0, 0, 0), 3: (ElemType.QUAD4, 4, 2, 1, 4, 0, 0, 0), 4: (ElemType.TETRA4, 4, 3, 1, 4, 0, 0, 0), 5: (ElemType.HEXA8, 8, 3, 1, 8, 0, 0, 0), 6: (ElemType.PRISM6, 6, 3, 1, 6, 0, 0, 0), 8: (ElemType.SEG3, 3, 1, 2, 2, 1, 0, 0), 9: (ElemType.TRI6, 6, 2, 2, 3, 3, 0, 0), 10: (ElemType.QUAD9, 9, 2, 2, 4, 4, 1, 0), 11: (ElemType.TETRA10, 10, 3, 2, 4, 6, 0, 0), 12: (ElemType.HEXA27, 27, 3, 2, 8, 12, 6, 1), 13: (ElemType.PRISM18, 18, 3, 2, 6, 9, 3, 0), 15: (ElemType.POINT, 1, 0, 0, 0, 0, 0, 0), 16: (ElemType.QUAD8, 8, 2, 2, 4, 4, 0, 0), 17: (ElemType.HEXA20, 20, 3, 2, 8, 12, 0, 0), 18: (ElemType.PRISM15, 15, 3, 2, 6, 9, 0, 0), 21: (ElemType.TRI10, 10, 2, 3, 3, 6, 1, 0), 23: (ElemType.TRI15, 15, 2, 4, 3, 9, 3, 0), 26: (ElemType.SEG4, 4, 1, 3, 2, 2, 0, 0), 27: (ElemType.SEG5, 5, 1, 4, 2, 3, 0, 0)}#

(ElemType, nPe, dim, order, Nvertex, Nedge, Nface, Nvolume)

Type:

gmshId

static Get_ElemInFos(gmshId)[source]#

return elemType, nPe, dim, order, Nvertex, Nedge, Nface, Nvolume

associated with the gmsh id.

Return type:

tuple[ElemType, int, int, int, int, int, int, int]

EasyFEA.fem.Inv(mat)[source]#

Computes inv(mat)

class EasyFEA.fem.LagrangeCondition(problemType, nodes, dofs, unknowns, dofsValues, lagrangeCoefs, description='')[source]#

Bases: BoundaryCondition

property lagrangeCoefs: ndarray[tuple[Any, ...], dtype[floating]]#

Lagrange coefficients.

class EasyFEA.fem.LinearForm(form)[source]#

Bases: _Form

Assemble(field)[source]#

Assemble de form with the field.

Parameters:

field (Field) – field

Returns:

the assembled (Nn * dof_n, 1) sparse vector.

Return type:

csr_matrix

Integrate_e(field)[source]#

Integrates de form with the field on elements.

Parameters:

field (Field) – field

Returns:

the integrated (Ne, …) numpy array

Return type:

np.ndarray

class EasyFEA.fem.MatrixType(value)[source]#

Bases: str, Enum

Order used for integration over elements, which determines the number of integration points.

static Get_types()[source]#
Return type:

list[MatrixType]

beam = 'beam'#

int_Ω ddNv • ddNv dΩ type

mass = 'mass'#

int_Ω N • N dΩ type

rigi = 'rigi'#

int_Ω dN • dN dΩ type

class EasyFEA.fem.Mesh(dict_groupElem, verbosity=False)[source]#

Bases: Observable

Mesh class that contains several _GroupElem instances.

Calc_regulation_projector(radius)[source]#

Returns the regulation projector matrix such that:

newU = proj * oldU

Parameters:

radius (float) – Regularization radius for the projection operation.

Returns:

Projection matrix of shape (Nn, Nn) that applies the regulation.

Return type:

sp.csr_matrix

Elements_Nodes(nodes, exclusively=True, neighborLayer=1)[source]#

Returns elements that exclusively or not use the specified nodes.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Elements_Tags(tags)[source]#

Returns elements associated with the tag.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Evaluate_dofsValues_at_coordinates(coordinates_n, dofsValues, elements=None)[source]#

Evaluates dofsValues with shape (Nn*dof_n, ) at the specified coordinates.

Parameters:
  • coordinates_n (_types.FloatArray) – coordinates that must be a (Nnodes, 3) array.

  • dofsValues (_types.FloatArray) – dofs values that must be a (Nn * dof_n) array.

  • elements (Optional[_types.IntArray], optional) – elements that may contain the specified coordinates to speed up evaluation, by default None

Returns:

The interpolated values as a (Nnodes, dof_n) array.

Return type:

_types.FloatArray

Get_B_e_pg(matrixType)[source]#

Get the matrix used to calculate deformations from displacements.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Warning

Use Kelvin Mandel Notation

[N1,x 0 … Nn,x 0

0 N1,y … 0 Nn,y

N1,y N1,x … N3,y N3,x]

(Ne, nPg, (3 or 6), nPe*dim)

Get_DiffusePart_e_pg(matrixType)[source]#

Get the part that builds the diffusion term (scalar).

DiffusePart_e_pg = k_e_pg * jacobian_e_pg * weight_pg * dN_e_pg’ @ A @ dN_e_pg

Returns (epij) -> jacobian_e_pg * weight_pg * dN_e_pg’

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_Gradient_e_pg(u, matrixType=MatrixType.rigi)[source]#

Returns the gradient of the discretized displacement field u as a matrix

Parameters:
  • u (_types.FloatArray) – discretized displacement field [ux1, uy1, uz1, …, uxN, uyN, uzN] of size Nn * dim

  • matrixType (MatrixType, optional) – matrix type, by default MatrixType.rigi

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Returns:

  • FeArray – grad(u) of shape (Ne, nPg, 3, 3)

  • dim = 1

  • ——-

  • dxux 0 0

  • 0 0 0

  • 0 0 0

  • dim = 2

  • ——-

  • dxux dyux 0

  • dxuy dyuy 0

  • 0 0 0

  • dim = 3

  • ——-

  • dxux dyux dzux

  • dxuy dyuy dzuy

  • dxuz dyuz dzuz

Get_N_pg(matrixType)[source]#

Evaluates shape functions in (ξ, η, ζ) coordinates.

[N1, … , Nn]

(nPg, 1, nPe)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_N_vector_pg(matrixType)[source]#

Returns shape functions matrix in (ξ, η, ζ) coordinates

[N1 0 … Nn 0

0 N1 … 0 Nn]

(nPg, dim, npe*dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_New_meshSize_n(error_e, coef=0.5)[source]#

Returns the scalar field (at nodes) used to refine the mesh.

meshSize = (coef - 1) / error_e.max() * error_e + 1

Parameters:
  • error_e (_types.FloatArray) – error evaluated on elements

  • coef (float, optional) – mesh size division ratio, by default 1/2

Returns:

meshSize_n, new mesh size at nodes (Nn)

Return type:

_types.FloatArray

Get_Node_Values(result_e)[source]#

Get node values from element values.

The value of a node is calculated by averaging the values of the surrounding elements.

Parameters:
  • mesh (Mesh) – mesh

  • result_e (_types.FloatArray) – element values (Ne, i)

Returns:

nodes values (Nn, i)

Return type:

_types.FloatArray

Get_Paired_Nodes(corners, plot=False)[source]#

Get the paired nodes used to construct periodic boundary conditions.

Parameters:
  • corners (_types.FloatArray) – Either nodes or nodes coordinates.

  • plot (bool, optional) – Set whether to plot the link between nodes; defaults to False.

Returns:

Paired nodes, a 2-column matrix storing the paired nodes (n, 2).

Return type:

_types.IntArray

Get_Quality(criteria='aspect', nodeValues=False)[source]#

Calculates mesh quality [0, 1] (bad, good).

Parameters:
  • criteria (str, optional) –

    criterion used, by default ‘aspect’

    • ”aspect”: hMin / hMax, ratio between minimum and maximum element length

    • ”angular”: angleMin / angleMax, ratio between the minimum and maximum angle of an element

    • ”gamma”: 2 rci/rcc, ratio between the radius of the inscribed circle and the circumscribed circle multiplied by 2. Useful for triangular elements.

    • ”jacobian”: jMax / jMin, ratio between the maximum jacobian and the minimum jacobian. Useful for higher-order elements.

  • nodeValues (bool, optional) – calculates values on nodes, by default False

Returns:

mesh quality

Return type:

_types.FloatArray

Get_ReactionPart_e_pg(matrixType)[source]#

Get the part that builds the reaction term (scalar).

ReactionPart_e_pg = r_e_pg * jacobian_e_pg * weight_pg * N_pg’ @ N_pg

Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’ @ N_pg

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_SourcePart_e_pg(matrixType)[source]#

Get the part that builds the source term (scalar).

SourcePart_e_pg = f_e_pg * jacobian_e_pg * weight_pg * N_pg’

Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_assembly_e(dof_n)[source]#

Returns assembly matrix for specified dof_n (Ne, nPe*dof_n)

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_connect_n_e()[source]#

Sparse matrix (Nn, Ne) of zeros and ones with ones when the node has the element such that:

values_n = connect_n_e * values_e

(Nn,1) = (Nn,Ne) * (Ne,1)

Return type:

csr_matrix

Get_dN_e_pg(matrixType)[source]#

Evaluates the first-order derivatives of shape functions in (x,y,z) coordinates.

[Ni,x … Nn,x

Ni,y … Nn,y]

(e, pg, dim, nPe)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_ddN_e_pg(matrixType)[source]#

Evaluates the first-order derivatives of shape functions in (x, y, z) coordinates.

[Ni,x … Nn,x

Ni,y … Nn,y

Ni,z … Nn,z]

(Ne, nPg, dim, nPe)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_jacobian_e_pg(matrixType, absoluteValues=True)[source]#

Returns the jacobians

variation in size (length, area or volume) between the reference element and the real element

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_leftDispPart(matrixType)[source]#

Get the left side of local displacement matrices.

Ku_e = jacobian_e_pg * weight_pg * B_e_pg’ @ c_e_pg @ B_e_pg

Returns (epij) -> jacobian_e_pg * weight_pg * B_e_pg’

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_list_groupElem(dim=None)[source]#

Returns the list of mesh element groups.

Parameters:

dim (int, optional) – The dimension of elements to retrieve, by default None (uses the main mesh dimension).

Returns:

A list of _GroupElem objects with the specified dimension.

Return type:

list[_GroupElem]

Get_meshSize(doMean=True)[source]#

Returns the mesh size of the mesh.

returns meshSize_e if doMean else meshSize_e_s

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_nPg(matrixType)[source]#

Returns integration points according to the matrix type.

Return type:

int

Get_normals(nodes=None, displacementMatrix=None)[source]#

Returns normal vectors and nodes belonging to the edge of the mesh.

returns normals, nodes.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]]

Get_weight_pg(matrixType)[source]#

Returns integration points according to the matrix type.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_weightedJacobian_e_pg(matrixType)[source]#

Returns jacobian_e_pg * weight_pg.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Locates_sol_e(sol, dof_n=None, asFeArray=False)[source]#

Locates solution on elements.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

property Ne: int#

number of elements in the mesh

property Nn: int#

number of nodes in the mesh

Nodes_Circle(circle, onlyOnEdge=True)[source]#

Returns the nodes in the circle.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Nodes_Conditions(func)[source]#

Returns nodes that meet the specified conditions.

Parameters:

func (function) –

Function using the x, y and z nodes coordinates and returning boolean values.

examples :

lambda x, y, z: (x < 40) & (x > 20) & (y<10)

lambda x, y, z: (x == 40) | (x == 50)

lambda x, y, z: x >= 0

Returns:

nodes that meet the specified conditions.

Return type:

_types.IntArray

Nodes_Cylinder(circle, direction=[0, 0, 1], onlyOnEdge=False)[source]#

Returns the nodes in the cylinder.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Nodes_Domain(domain)[source]#

Returns nodes in the domain.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Nodes_Line(line)[source]#

Returns the nodes on the line.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Nodes_Point(point)[source]#

Returns nodes on the point.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Nodes_Points(points)[source]#

Returns nodes on points.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Nodes_Tags(tags)[source]#

Returns nodes associated with the tags.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Rotate(theta, center=(0, 0, 0), direction=(0, 0, 1))[source]#

Rotates the mesh coordinates around an axis.

Parameters:
  • theta (float) – rotation angle [deg]

  • center (tuple, optional) – rotation center, by default (0,0,0)

  • direction (tuple, optional) – rotation direction, by default (0,0,1)

Return type:

None

Set_Tag(nodes, tag)[source]#

Set a tag on the nodes and elements belonging to each group of elements in the mesh.

Symmetry(point=(0, 0, 0), n=(1, 0, 0))[source]#

Symmetrizes the mesh coordinates with respect to a specified plane.

Parameters:
  • point (tuple, optional) – a point belonging to the plane, by default (0,0,0)

  • n (tuple, optional) – normal to the plane, by default (1,0,0)

Return type:

None

Translate(dx=0.0, dy=0.0, dz=0.0)[source]#

Translates the mesh coordinates.

Return type:

None

property area: float#

total area of the mesh.

property assembly_e: ndarray[tuple[Any, ...], dtype[IntType]]#

assembly matrix (Ne, nPe*dim)

Used to position the rigi matrix in the global matrix.

property center: ndarray[tuple[Any, ...], dtype[floating]]#

center of mass / barycenter / inertia center

property columnsScalar_e: ndarray[tuple[Any, ...], dtype[IntType]]#

columns to fill the assembly matrix in scalar form (damage or thermal problems)

property columnsVector_e: ndarray[tuple[Any, ...], dtype[IntType]]#

columns to fill the assembly matrix in vector (e.g elastic problem)

property connect: ndarray[tuple[Any, ...], dtype[IntType]]#

connectivity matrix (Ne, nPe)

property coord: ndarray[tuple[Any, ...], dtype[floating]]#

global nodes coordinates matrix (Nn, 3)

Contains all nodes coordinates

copy()[source]#
property dict_groupElem: dict[ElemType, _GroupElem]#

dictionary containing all the element groups in the mesh

property dim#

mesh dimension

property elemType: ElemType#

mesh element type

property groupElem: _GroupElem#

main group element

property inDim#

dimension in which the mesh lies.

A 2D mesh can be oriented in a 3D space.

property length: float#

total length of the mesh.

property nPe: int#

nodes per element

property nodes: ndarray[tuple[Any, ...], dtype[IntType]]#

mesh nodes

property orphanNodes: list[int]#

nodes not connected to any mesh elements

property rowsScalar_e: ndarray[tuple[Any, ...], dtype[IntType]]#

rows to fill the assembly matrix in scalar form (damage or thermal problems)

property rowsVector_e: ndarray[tuple[Any, ...], dtype[IntType]]#

rows to fill the assembly matrix in vector (e.g elastic problem)

property verbosity: bool#

the mesh can write in the terminal

property volume: float#

total volume of the mesh.

EasyFEA.fem.Mesh_Optim(DoMesh, folder, criteria='aspect', quality=0.8, ratio=0.7, iterMax=20, coef=0.5)[source]#

Optimize the mesh using the given criterion.

Parameters:
  • DoMesh (Callable[[str], Mesh]) –

    Function that constructs the mesh and takes a .pos file as argument for mesh optimization.

    The function must return a Mesh.

  • folder (str) – Folder in which .pos files are created and then deleted.

  • criteria (str, optional) –

    criterion used, by default ‘aspect’

    • ”aspect”: hMin / hMax, ratio between minimum and maximum element length

    • ”angular”: angleMin / angleMax, ratio between the minimum and maximum angle of an element

    • ”gamma”: 2 rci/rcc, ratio between the radius of the inscribed circle and the circumscribed circle multiplied by 2. Useful for triangular elements.

    • ”jacobian”: jMax / jMin, ratio between the maximum jacobian and the minimum jacobian. Useful for higher-order elements.

  • quality (float, optional) – quality target, by default .8

  • ratio (float, optional) – target ratio of mesh elements that must respect the specified quality, by default 0.7 (must be in [0,1])

  • iterMax (int, optional) – Maximum number of iterations, by default 20

  • coef (float, optional) – mesh size division ratio, by default 1/2

Returns:

optimized mesh size and ratio

Return type:

tuple[Mesh, float]

class EasyFEA.fem.Mesher(openGmsh=False, gmshVerbosity=False, verbosity=False)[source]#

Bases: object

Mesher class used to construct and generate the mesh via gmsh.

Create_posFile(coord, values, folder, filename='data')[source]#

Creates of a .pos file that can be used to refine a mesh in a zone.

Parameters:
  • coord (_types.FloatArray) – coordinates of values

  • values (_types.FloatArray) – scalar nodes values

  • folder (str) – save folder

  • filename (str, optional) – .pos file name, by default “data”.

Returns:

Returns the path to the created .pos file

Return type:

str

static Get_Entities(points=[], lines=[], surfaces=[], volumes=[])[source]#

Get entities from from points, lines, surfaces and volumes tags

Return type:

list[tuple[int, int]]

Mesh_2D(contour, inclusions=[], elemType=ElemType.TRI3, cracks=[], refineGeoms=[], isOrganised=False, additionalSurfaces=[], additionalLines=[], additionalPoints=[], folder='')[source]#

Creates a 2D mesh from a contour and inclusions that must form a closed plane surface.

Parameters:
  • contour (Domain | Circle | Points | Contour) – geom object

  • inclusions (list[Domain | Circle | Points | Contour], optional) – list of hollow and filled geom objects inside the domain

  • elemType (ElemType, optional) – element type, by default “TRI3” [“TRI3”, “TRI6”, “TRI10”, “TRI15”, “QUAD4”, “QUAD8”, “QUAD9”]

  • cracks (list[Line | Points | Contour | CircleArc]) – list of geom object used to create open or closed cracks

  • refineGeoms (list[Domain|Circle|str], optional) – list of geom object for mesh refinement, by default []

  • isOrganised (bool, optional) – mesh is organized, by default False

  • additionalSurfaces (list[Domain | Circle | Points | Contour]) – additional surfaces that will be added to or removed from the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). Tip: if the mesh is not well generated, you can also give the inclusions.

  • additionalLines (list[Union[Line,CircleArc]]) – additional lines that will be added to the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). WARNING: lines must be within the domain.

  • additionalPoints (list[Point]) – additional points that will be added to the surfaces created by the contour and the inclusions. WARNING: points must be within the domain.

  • folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Beams(beams, elemType=ElemType.SEG2, additionalPoints=[], folder='')[source]#

Creates a beam mesh.

Parameters:
  • beams (list[_Beam]) – list of Beams

  • elemType (ElemType, optional) – element type, by default “SEG2” [“SEG2”, “SEG3”, “SEG4”]

  • folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh

  • additionalPoints (list[Point]) – additional points that will be added to the mesh. WARNING: points must be within the domain.

Returns:

Created mesh

Return type:

Mesh

Mesh_Extrude(contour, inclusions=[], extrude=(0, 0, 1), layers=[], elemType=ElemType.TETRA4, cracks=[], refineGeoms=[], isOrganised=False, additionalSurfaces=[], additionalLines=[], additionalPoints=[], folder='')[source]#

Creates a 3D mesh by extruding a surface constructed from a contour and inclusions.

Parameters:
  • contour (Domain | Circle | Points | Contour) – geom object

  • inclusions (list[Domain | Circle | Points | Contour], optional) – list of hollow and filled geom objects inside the domain

  • extrude (Coords, optional) – extrusion vector, by default [0,0,1]

  • layers (list[int], optional) – layers in the extrusion, by default []

  • elemType (ElemType, optional) – element type, by default “TETRA4” [“TETRA4”, “TETRA10”, “HEXA8”, “HEXA20”, “HEXA27”, “PRISM6”, “PRISM15”, “PRISM18”]

  • cracks (list[Line | Points | Contour | CircleArc]) – list of geom object used to create open or closed cracks

  • refineGeoms (list[Domain|Circle|str], optional) – list of geom object for mesh refinement, by default []

  • isOrganised (bool, optional) – mesh is organized, by default False

  • additionalSurfaces (list[Domain | Circle | Points | Contour]) – additional surfaces that will be added to or removed from the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). Tip: if the mesh is not well generated, you can also give the inclusions.

  • additionalLines (list[Union[Line,CircleArc]]) – additional lines that will be added to the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). WARNING: lines must be within the domain.

  • additionalPoints (list[Point]) – additional points that will be added to the surfaces created by the contour and the inclusions. WARNING: points must be within the domain.

  • folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Import_mesh(mesh, setPhysicalGroups=False, coef=1.0)[source]#

Creates the mesh from an .msh file.

Parameters:
  • mesh (str) – .msh file

  • setPhysicalGroups (bool, optional) – creates physical groups, by default False

  • coef (float, optional) – coef applied to the node coordinates, by default 1.0

Returns:

Created mesh

Return type:

Mesh

Mesh_Import_part(file, dim, meshSize=0.0, elemType=None, refineGeoms=[None], folder='')[source]#

Creates the mesh from .stp or .igs files.

You can only use triangles or tetrahedrons.

Parameters:
  • file (str) –

    .stp or .igs files.

    Note that for igs files, entities cannot be recovered.

  • meshSize (float, optional) – mesh size, by default 0.0

  • elemType (ElemType, optional) – element type, by default “TRI3” or “TETRA4” depending on dim.

  • refineGeoms (list[Domain|Circle|str]) – geometric objects to refine the mesh

  • folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Revolve(contour, inclusions=[], axis=<EasyFEA.geoms._line.Line object>, angle=360, layers=[30], elemType=ElemType.TETRA4, cracks=[], refineGeoms=[], isOrganised=False, additionalSurfaces=[], additionalLines=[], additionalPoints=[], folder='')[source]#

Creates a 3D mesh by rotating a surface along an axis.

Parameters:
  • contour (Domain | Circle | Points | Contour) – geometry that builds the contour

  • inclusions (list[Domain | Circle | Points | Contour], optional) – list of hollow and filled geom objects inside the domain

  • axis (Line, optional) – revolution axis, by default Line(Point(), Point(0,1))

  • angle (float|int, optional) – revolution angle in [deg], by default 360

  • layers (list[int], optional) – layers in extrusion, by default [30]

  • elemType (ElemType, optional) – element type, by default “TETRA4” [“TETRA4”, “TETRA10”, “HEXA8”, “HEXA20”, “HEXA27”, “PRISM6”, “PRISM15”, “PRISM18”]

  • cracks (list[Line | Points | Contour | CircleArc]) – list of geom object used to create open or closed cracks

  • refineGeoms (list[Domain|Circle|str], optional) – list of geom object for mesh refinement, by default []

  • isOrganised (bool, optional) – mesh is organized, by default False

  • additionalSurfaces (list[Domain | Circle | Points | Contour]) – additional surfaces that will be added to or removed from the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). Tip: if the mesh is not well generated, you can also give the inclusions.

  • additionalLines (list[Union[Line,CircleArc]]) – additional lines that will be added to the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). WARNING: lines must be within the domain.

  • additionalPoints (list[Point]) – additional points that will be added to the surfaces created by the contour and the inclusions. WARNING: points must be within the domain.

  • folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Save_Simu(simu, results=[], details=False, edgeColor='black', plotMesh=True, showAxes=False, folder='')[source]#

Save the simulation in gmsh.pos format using gmsh.view

Parameters:
  • simu (_Simu) – simulation

  • results (list[str], optional) – list of result you want to plot, by default []

  • details (bool, optional) – get default result values with details or not see simu.Results_nodesField_elementsField(details), by default False

  • edgeColor (str, optional) – color used to plot the edges, by default ‘black’

  • plotMesh (bool, optional) – plot the mesh, by default True

  • showAxes (bool, optional) – show the axes, by default False

  • folder (str, optional) – folder used to save .pos file, by default “”

Return type:

None

Set_meshSize(meshSize)[source]#

Sets the mesh size

Return type:

None

EasyFEA.fem.Norm(array, **kwargs)[source]#

np.linalg.norm() wrapper.

see https://numpy.org/doc/stable/reference/generated/numpy.linalg.norm.html

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

EasyFEA.fem.Sym_Grad(u)[source]#
Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

EasyFEA.fem.TensorProd(A, B, symmetric=False, ndim=None)[source]#

Computes tensor product.

Parameters:
  • A (FeArray.FeArrayALike) – array A

  • B (FeArray.FeArrayALike) – array B

  • symmetric (bool, optional) – do symmetric product, by default False

  • ndim (int, optional) – ndim=1 -> vect or ndim=2 -> matrix, by default None

Returns:

the calculated tensor product

Return type:

FeArray.FeArrayALike

EasyFEA.fem.Trace(mat)[source]#

Computes trace(mat)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

EasyFEA.fem.Transpose(mat)[source]#

Computes transpose(mat)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

class EasyFEA.fem._GroupElem(gmshId, connect, coordGlob)[source]#

Bases: ABC

The _GroupElem base class, from which all element types inherit.

Get_B_e_pg(matrixType)[source]#

Get the matrix used to calculate deformations from displacements.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Warning

Use Kelvin Mandel Notation

[N1,x 0 … Nn,x 0

0 N1,y … 0 Nn,y

N1,y N1,x … N3,y N3,x]

(Ne, nPg, (3 or 6), nPe*dim)

Get_DiffusePart_e_pg(matrixType)[source]#

Get the part that builds the diffusion term (scalar).

DiffusePart_e_pg = k_e_pg * jacobian_e_pg * weight_pg * dN_e_pg’ @ A @ dN_e_pg

Returns (epij) -> jacobian_e_pg * weight_pg * dN_e_pg’

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_Elements_Nodes(nodes, exclusively=True)[source]#

Returns elements that exclusively or not use the specified nodes.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_Elements_Tag(tag)[source]#

Returns elements associated with the tag.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_EulerBernoulli_B_e_pg(beamStructure)[source]#

Get Euler-Bernoulli beam shape functions derivatives

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_EulerBernoulli_N_e_pg()[source]#

Evaluates Euler-Bernoulli beam shape functions in (x, y, z) coordinates.

[phi_i psi_i … phi_n psi_n]

(Ne, nPg, 1, nPe*2)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_EulerBernoulli_N_e_pg_for_beam(beamStructure)[source]#

Euler-Bernoulli beam shape functions.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_EulerBernoulli_N_pg()[source]#

Evaluates Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.

[phi_i psi_i … phi_n psi_n]

(nPg, nPe*2)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_EulerBernoulli_dN_e_pg()[source]#

Evaluates the first-order derivatives of Euler-Bernoulli beam shape functions in (x, y, z) coordinates.

[phi_i,x psi_i,x … phi_n,x psi_n,x]

(Ne, nPg, 1, nPe*2)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_EulerBernoulli_dN_pg()[source]#

Evaluates Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPg, nPe*2)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_EulerBernoulli_ddN_e_pg()[source]#

Evaluates the second-order derivatives of Euler-Bernoulli beam shape functions in (x, y, z) coordinates.

[phi_i,xx psi_i,xx … phi_n,xx psi_n,xx]

(Ne, nPg, 1, nPe*2)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_EulerBernoulli_ddN_pg()[source]#

Evaluates Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ x psi_n,ξ]

(nPg, nPe*2)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_F_e_pg(matrixType)[source]#

Returns the transposed Jacobian matrix.

This matrix describes the transformation of the (ξ, η, ζ) axes from the reference element to the (x, y, z) coordinate system of the actual element.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_GaussCoordinates_e_pg(matrixType, elements=array([], dtype=float64))[source]#

Returns integration point coordinates for each element (Ne, nPg, 3) in the (x, y, z) coordinates.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_Gradient_e_pg(u, matrixType=MatrixType.rigi)[source]#

Returns the gradient of the discretized displacement field u as a matrix

Parameters:
  • u (_types.FloatArray) – discretized displacement field [ux1, uy1, uz1, …, uxN, uyN, uzN] of size Nn * dim

  • matrixType (MatrixType, optional) – matrix type, by default MatrixType.rigi

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Returns:

  • FeArray – grad(u) of shape (Ne, nPg, 3, 3)

  • dim = 1

  • ——-

  • dxux 0 0

  • 0 0 0

  • 0 0 0

  • dim = 2

  • ——-

  • dxux dyux 0

  • dxuy dyuy 0

  • 0 0 0

  • dim = 3

  • ——-

  • dxux dyux dzux

  • dxuy dyuy dzuy

  • dxuz dyuz dzuz

abstractmethod Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_Mapping(coordinates_n, elements_e=None, needCoordinates=False)[source]#

Locates coordinates within elements.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]], ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]], ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]], Optional[ndarray[tuple[Any, ...], dtype[floating]]]]

Returns:

  • detectedNodes : The nodes that have been identified within the detected elements with shape=(Nn,).

  • detectedElements_e : The elements in which the nodes have been detected with shape=(Ne,).

  • connect_e_nThe connectivity matrix that includes the nodes identified in each element with shape=(Ne, ?).

    The “?” indicates that the array may have varying dimensions along axis=1.

  • coordInElem_nThe coordinates of the identified nodes, expressed in the reference element’s (ξ, η, ζ) coordinate system.

    This is applicable only if needCoordinates is set to True.

Get_N_pg(matrixType)[source]#

Evaluates shape functions in (ξ, η, ζ) coordinates.

[N1, … , Nn]

(nPg, 1, nPe)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_N_pg_rep(matrixType, repeat=1)[source]#

Repeats shape functions in the (ξ, η, ζ) coordinates.

Parameters:
  • matrixType (MatrixType) – matrix type

  • repeat (int, optional) – number of repetitions, by default 1

  • Returns

  • -------

  • (nPg (• Scalar shape functions) –

    [Ni 0 … Nn 0

    0 Ni … 0 Nn]

  • rep=2

    [Ni 0 … Nn 0

    0 Ni … 0 Nn]

  • rep=2*dim)

    [Ni 0 … Nn 0

    0 Ni … 0 Nn]

  • (nPg – [Ni … Nn]

  • rep=1 – [Ni … Nn]

  • nPe) – [Ni … Nn]

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_Nodes_Circle(circle, onlyOnEdge=False)[source]#

Returns nodes in the circle.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_Nodes_Conditions(func)[source]#

Returns nodes that meet the specified conditions.

Parameters:

func (Callable) –

Function using x, y and z nodes coordinates and returning boolean values.

examples :

lambda x, y, z: (x < 40) & (x > 20) & (y<10)

lambda x, y, z: (x == 40) | (x == 50)

lambda x, y, z: x >= 0

Returns:

nodes that meet conditions

Return type:

_types.IntArray

Get_Nodes_Cylinder(circle, direction=[0, 0, 1], onlyOnEdge=False)[source]#

Returns nodes in the cylinder.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_Nodes_Domain(domain)[source]#

Returns nodes in the domain.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_Nodes_Line(line)[source]#

Returns nodes on the line.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_Nodes_Point(point)[source]#

Returns nodes on the point.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_Nodes_Tag(tag)[source]#

Returns node associated with the tag.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_ReactionPart_e_pg(matrixType)[source]#

Get the part that builds the reaction term (scalar).

ReactionPart_e_pg = r_e_pg * jacobian_e_pg * weight_pg * N_pg’ @ N_pg

Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’ @ N_pg

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_SourcePart_e_pg(matrixType)[source]#

Get the part that builds the source term (scalar).

SourcePart_e_pg = f_e_pg * jacobian_e_pg * weight_pg * N_pg’

Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_assembly_e(dof_n)[source]#

Get the assembly matrix for the specified dof_n (Ne, nPe*dof_n)

Parameters:

dof_n (int) – degree of freedom per node

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_columnsVector_e(dof_n)[source]#

Returns columns to fill the vector assembly matrix

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_connect_n_e()[source]#

Sparse matrix (Nn, Ne) of zeros and ones with ones when the node has the element such that: values_n = connect_n_e * values_e

(Nn,1) = (Nn,Ne) * (Ne,1)

Return type:

csr_matrix

Get_dN_e_pg(matrixType)[source]#

Evaluates the first-order derivatives of shape functions in (x, y, z) coordinates.

[Ni,x … Nn,x

Ni,y … Nn,y

Ni,z … Nn,z]

(Ne, nPg, dim, nPe)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_dN_pg(matrixType)[source]#

Evaluates shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ … Nn,ξ

Ni,η … Nn,η

Ni,ζ … Nn,ζ

(nPg, dim, nPe)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_ddN_e_pg(matrixType)[source]#

Evaluates the second-order derivatives of shape functions in (x, y, z) coordinates.

[Ni,x2 … Nn,x2

Ni,y2 … Nn,y2

Ni,z2 … Nn,z2]

(Ne, nPg, dim, nPe)

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_ddN_pg(matrixType)[source]#

Evaluates shape functions second derivatives in the (ξ, η, ζ) coordinates.

[Ni,ξ2 … Nn,ξ2

Ni,η2 … Nn,η2

Ni,ζ2 … Nn,ζ2]

(nPg, dim, nPe)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_dddN_pg(matrixType)[source]#

Evaluates shape functions third derivatives in the (ξ, η, ζ) coordinates.

[Ni,ξ3 … Nn,ξ3

Ni,η3 … Nn,η3

Ni,ζ3 … Nn,ζ3]

(nPg, dim, nPe)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_ddddN_pg(matrixType)[source]#

Evaluates shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

[Ni,ξ4 … Nn,ξ4

Ni,η4 … Nn,η4

Ni,ζ4 … Nn,ζ4]

(pg, dim, nPe)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_gauss(matrixType)[source]#

Returns integration points according to the matrix type.

Return type:

Gauss

Get_invF_e_pg(matrixType)[source]#

Returns the inverse of the transposed Jacobian matrix.

Used to obtain the derivative of the dN_e_pg shape functions in the actual element dN_e_pg = invF_e_pg • dN_pg

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_jacobian_e_pg(matrixType, absoluteValues=True)[source]#

Returns the jacobians.

variation in size (length, area or volume) between the reference element and the actual element

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_leftDispPart(matrixType)[source]#

Get the left side of local displacement matrices.

Ku_e = jacobian_e_pg * weight_pg * B_e_pg’ @ c_e_pg @ B_e_pg

Returns (epij) -> jacobian_e_pg * weight_pg * B_e_pg’

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Get_normals_e_pg(matrixType, displacementMatrix=None)[source]#

Returns the normals for each elements and gauss points (Ne, nPg, 3).

Return type:

FeArray

Get_pointsInElem(coordinates_n, elem)[source]#

Returns the indexes of the coordinates contained in the element.

Parameters:
  • coordinates_n (_types.FloatArray) – coordinates (n, 3)

  • elem (int) – element

Returns:

indexes of coordinates contained in element

Return type:

_types.IntArray

Get_rowsVector_e(dof_n)[source]#

Returns rows to fill the assembly matrix in vector (e.g. elastic problem)

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Get_weight_pg(matrixType)[source]#

Returns integration point weights according to the matrix type.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Get_weightedJacobian_e_pg(matrixType)[source]#

Returns the jacobian_e_pg * weight_pg.

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

Integrate_e(func=<function _GroupElem.<lambda>>, matrixType=MatrixType.mass)[source]#

Integrates the function over elements.

Parameters:
  • func (lambda) –

    function that uses the x,y,z coordinates of the element’s integration points

    Examples:

    lambda x,y,z: 1 -> that will just integrate the element lambda x,y,z: x lambda x,y,z: x + y

    lambda x,y,z: z**2

  • matrixType (MatrixType, optional) – matrix type, by default MatrixType.mass

Returns:

integrated values on elements

Return type:

_types.FloatArray

Locates_sol_e(sol, dof_n=None, asFeArray=False)[source]#

Locates sol on elements

Return type:

Union[FeArray, ndarray[tuple[Any, ...], dtype[Any]]]

property Ne: int#

number of elements

property Nedge: int#

number of edge nodes per element

property Nface: int#

number of face nodes per element

property Nn: int#

number of nodes used by the element group

property Nvertex: int#

number of vertex nodes per element

property Nvolume: int#

number of volume nodes per element

Set_Tag(nodes, tag)[source]#

Set a tag on the nodes and elements belonging to the group of elements.

property area: float#

area covered by elements

property area_e: ndarray[tuple[Any, ...], dtype[floating]]#

area covered by each element

property assembly_e: ndarray[tuple[Any, ...], dtype[IntType]]#

assembly matrix (Ne, nPe*dim)

property center: ndarray[tuple[Any, ...], dtype[floating]]#

center of mass / barycenter / inertia center

property connect: ndarray[tuple[Any, ...], dtype[IntType]]#

connectivity matrix (Ne, nPe)

property coord: ndarray[tuple[Any, ...], dtype[floating]]#

this matrix contains the element group coordinates (Nn, 3)

property coordGlob: ndarray[tuple[Any, ...], dtype[floating]]#

this matrix contains all the mesh coordinates (mesh.Nn, 3)

property dim: int#

element dimension

property edges: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element edges (for FEM purposes).

property elemType: ElemType#

element type

property elementTags: list[str]#

returns element tags.

property elements: ndarray[tuple[Any, ...], dtype[IntType]]#
abstract property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property gmshId: int#

gmsh Id

property inDim: int#

dimension in which the elements are located

property length: float#

length covered by elements

property length_e: ndarray[tuple[Any, ...], dtype[floating]]#

length covered by each element

property nPe: int#

nodes per element

property nodeTags: list[str]#

Returns node tags.

property nodes: ndarray[tuple[Any, ...], dtype[IntType]]#

nodes used by the element group. Node ‘n’ is on line ‘n’ in coordGlob

property order: int#

element order

abstract property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

abstract property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property topology: str#

element topology

abstract property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

property volume: float#

volume covered by elements

property volume_e: ndarray[tuple[Any, ...], dtype[floating]]#

volume covered by each element

class EasyFEA.fem.elems.HEXA20(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.HEXA27(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.HEXA8(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.POINT(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.PRISM15(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.PRISM18(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.PRISM6(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.QUAD4(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.QUAD8(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.QUAD9(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.SEG2(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_EulerBernoulli_N()[source]#

Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_dN()[source]#

Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_ddN()[source]#

Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 2)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.SEG3(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_EulerBernoulli_N()[source]#

Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_dN()[source]#

Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_ddN()[source]#

Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 2)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.SEG4(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_EulerBernoulli_N()[source]#

Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_dN()[source]#

Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_ddN()[source]#

Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 2)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.SEG5(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_EulerBernoulli_N()[source]#

Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_dN()[source]#

Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 1)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_EulerBernoulli_ddN()[source]#

Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.

[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]

(nPe*2, 2)

Return type:

ndarray[tuple[Any, ...], dtype[Any]]

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.TETRA10(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.TETRA4(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property segments: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to construct segments (for display purposes).

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.TRI10(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.TRI15(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.TRI3(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function

class EasyFEA.fem.elems.TRI6(gmshId, connect, coordGlob)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array

_N()[source]#

Shape functions in (ξ, η, ζ) coordinates.

N1

Nn

(nPe, 1)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_abc_impl = <_abc._abc_data object>#
_dN()[source]#

Shape functions first derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddN()[source]#

Shape functions second derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ2 Ni,η2 Ni,ζ2

Nn,ξ2 Nn,η2 Nn,ζ2

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_dddN()[source]#

Shape functions third derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ3 Ni,η3 Ni,ζ3

Nn,ξ3 Nn,η3 Nn,ζ3

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

_ddddN()[source]#

Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.

Ni,ξ4 Ni,η4 Ni,ζ4

Nn,ξ4 Nn,η4 Nn,ζ4

(nPe, dim)

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

property faces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the element faces (for FEM purposes).

property origin: list[int]#

reference element origin coordinates

property surfaces: ndarray[tuple[Any, ...], dtype[IntType]]#

array of indices used to form the contour of the surfaces that make up the element (for display purposes).

Warning

When adding new 3D elements, ensure that the resulting surface normals point inward the element.

property triangles: list[int]#

list of index used to form the triangles of an element that will be used for the 2D trisurf function