FEM#

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

In the simulation workflow, FEM is the second step: geometry objects from Geoms are passed to Mesher to produce the Mesh.

See also


What is a mesh in EasyFEA?#

A Mesh object in EasyFEA represents a collection of elements 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 module and were defined in ShapeFunctions. The Gauss point quadratures are implemented in the Gauss class.


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:

FEM API#

class EasyFEA.FEM.BiLinearForm(
form: Callable[[...], FeArray | ndarray[tuple[Any, ...], dtype[Any]]],
)[source]#

Bases: _Form

Bilinear form.

Assemble(field) csr_matrix[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: str,
nodes: ndarray[tuple[Any, ...], dtype[int64]],
dofs: ndarray[tuple[Any, ...], dtype[int64]],
unknowns: list[str],
dofsValues: ndarray[tuple[Any, ...], dtype[floating]],
description: str,
)[source]#

Bases: object

Boundary condition.

static Get_dofs(
problemType: str,
list_Bc_Condition: list[BoundaryCondition],
) ndarray[tuple[Any, ...], dtype[int64]][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: list[str],
nodes: ndarray[tuple[Any, ...], dtype[int64]],
unknowns: list[str],
) ndarray[tuple[Any, ...], dtype[int64]][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: str,
list_Bc_Condition: list[BoundaryCondition],
) int[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: str,
list_Bc_Condition: list[BoundaryCondition],
) ndarray[tuple[Any, ...], dtype[floating]][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[int64]][source]#

degrees of freedom associated with the nodes and unknowns

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

values applied

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

nodes on which the condition is applied

property problemType: str[source]#

type of problem

property unknowns: list[str][source]#

dofs unknowns

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

Bases: str, Enum

Implemented Lagrange isoparametric element types.

static Get_1D() list[ElemType][source]#

Returns 1D element types.

static Get_2D() list[ElemType][source]#

Returns 2D element types.

static Get_3D() list[ElemType][source]#

Returns 3D element types.

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[source]#
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.

static _ddot_subscript(ndim1: int, ndim2: int) str[source]#

Build and cache the einsum subscript for ddot(ndim1, ndim2).

static _dot_subscript(ndim1: int, ndim2: int) str[source]#

Build and cache the einsum subscript for dot(ndim1, ndim2).

static asfearray(
array,
broadcastFeArrays=False,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
static ones(
*shape,
dtype=None,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
static zeros(
*shape,
dtype=None,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
_asfearrays(
*,
broadcastFeArrays=False,
) list[FeArray | ndarray[tuple[Any, ...], dtype[Any]]][source]#
_assemble(
*arrays,
value: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
)[source]#
_get_idx(
*arrays,
) list[ndarray[tuple[Any, ...], dtype[Any]]][source]#
ddot(
other,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
dot(other, /, out=None)[source]#

Refer to numpy.dot() for full documentation.

See also

numpy.dot

equivalent function

max(
*args,
**kwargs,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

np.max() wrapper.

min(
*args,
**kwargs,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

np.min() wrapper.

sum(
*args,
**kwargs,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

np.sum() wrapper.

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

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

property _idx: str[source]#

einsum indicator (e.g “”, “i”, “ij”) used in np.einsum() function.

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

property _ndim: int[source]#

finite element ndim

property _shape: tuple[source]#

finite element shape

property _type: str[source]#
class EasyFEA.FEM.Field(
groupElem: _GroupElem,
dof_n: int,
matrixType: MatrixType = MatrixType.mass,
)[source]#

Bases: object

Field class.

Evaluate_e(
function: Callable[[Field], FeArray],
dofsValues: ndarray,
returnMeanValues: bool = True,
) ndarray[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.

Evaluate_n(
function: Callable[[Field], FeArray],
dofsValues: ndarray,
) ndarray[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.

Get_coords(concatenate=False)[source]#

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

Interpolate(dofsValues: ndarray) FeArray[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

_Get_current_active_dof() int[source]#

Returns current active dof.

_Get_current_active_node() int[source]#

Returns current active node.

_Get_dofsValues() ndarray[tuple[Any, ...], dtype[floating]][source]#
_Set_current_active_dof(dof: int)[source]#

Sets current active node.

_Set_current_active_node(node: int)[source]#

Sets current active node.

_Set_dofsValues(
values: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#
copy() Field[source]#
ddot(other)[source]#
dot(other)[source]#
property dof_n: int[source]#

degrees of freedom per node.

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

Returns the gradient of the field.

property groupElem: _GroupElem[source]#

Group of elements.

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

Bases: object

Gauss quadrature.

static Gauss_factory(
elemType: ElemType,
matrixType: MatrixType,
) tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#

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

static _Gauss_factory_nPg(
elemType: ElemType,
nPg,
) tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#

Returns the integration points and weights based on the element type and the number of Gauss points (nPg).

static _Hexahedron(
nPg: int,
) tuple[Collection[int | float], Collection[int | float], Collection[int | float], Collection[int | float]][source]#

available [8, 27]

order = [3, 5]

static _Prism(
nPg: int,
) tuple[Collection[int | float], Collection[int | float], Collection[int | float], Collection[int | float]][source]#

available [6, 8, 21]

order X = [3, 3, 5]

order Y & Z = [2, 3, 5]

static _Quadrangle(
nPg: int,
) tuple[Collection[int | float], Collection[int | float], Collection[int | float]][source]#

available [4, 9]

order = [1, 2]

static _Tetrahedron(
nPg: int,
) tuple[Collection[int | float], Collection[int | float], Collection[int | float], Collection[int | float]][source]#

available [1, 4, 5, 15]

order = [1, 2, 3, 5]

static _Triangle(
nPg: int,
) tuple[Collection[int | float], Collection[int | float], Collection[int | float]][source]#

available [1, 3, 6, 7, 12]

order = [1, 2, 3, 4, 5]

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

integration point coordinates

property nPg: int[source]#

number of integration points

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

integration point weights

class EasyFEA.FEM.GroupElemFactory[source]#

Bases: object

static Create(
elemType: ElemType,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
) _GroupElem[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

static Get_ElemInFos(
gmshId: int,
) tuple[ElemType, int, int, int, int, int, int, int][source]#

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

associated with the gmsh id.

static _Create(
gmshId: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
) _GroupElem[source]#

Creates an element group.

Parameters:
  • gmshId (int) – id gmsh

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

  • coordGlob (_types.FloatArray) – nodes coordinates

Returns:

the element group

Return type:

GroupeElem

static _Get_2d_element_types(
elemType: ElemType,
) list[ElemType][source]#

Returns 2d element types associated with the elementType.

Parameters:

elemType (ElemType) – element type

Returns:

the element types

Return type:

list[ElemType]

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

class EasyFEA.FEM.LagrangeCondition(
problemType: str,
nodes: ndarray[tuple[Any, ...], dtype[int64]],
dofs: ndarray[tuple[Any, ...], dtype[int64]],
unknowns: list[str],
dofsValues: ndarray[tuple[Any, ...], dtype[floating]],
lagrangeCoefs: ndarray[tuple[Any, ...], dtype[floating]],
description: str = '',
)[source]#

Bases: BoundaryCondition

Lagrange boundary condition

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

Lagrange coefficients.

class EasyFEA.FEM.LinearForm(
form: Callable[[...], FeArray | ndarray[tuple[Any, ...], dtype[Any]]],
)[source]#

Bases: _Form

Linear form.

Assemble(field) csr_matrix[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() list[MatrixType][source]#
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: dict[ElemType, _GroupElem], verbosity=False)[source]#

Bases: Observable

Mesh class that contains several _GroupElem instances.

Calc_regulation_projector(radius: float) csr_matrix[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: ndarray[tuple[Any, ...], dtype[int64]],
exclusively=True,
neighborLayer: int = 1,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns elements that exclusively or not use the specified nodes.

Elements_Tags(
*tags: str,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns elements associated with the tag.

Evaluate_dofsValues_at_coordinates(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
dofsValues: ndarray[tuple[Any, ...], dtype[floating]],
elements: ndarray[tuple[Any, ...], dtype[int64]] | None = None,
) ndarray[tuple[Any, ...], dtype[floating]][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: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Get the matrix used to calculate deformations from displacements.

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: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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’

Get_Gradient_e_pg(
u: ndarray[tuple[Any, ...], dtype[floating]],
matrixType=MatrixType.rigi,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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

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: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[N1, … , Nn]

(nPg, 1, nPe)

Get_N_vector_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[N1 0 … Nn 0

0 N1 … 0 Nn]

(nPg, dim, npe*dim)

Get_New_meshSize_n(
error_e: ndarray[tuple[Any, ...], dtype[floating]],
coef: float = 0.5,
) ndarray[tuple[Any, ...], dtype[floating]][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: ndarray[tuple[Any, ...], dtype[floating]],
) ndarray[tuple[Any, ...], dtype[floating]][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: ndarray[tuple[Any, ...], dtype[floating]],
plot=False,
) ndarray[tuple[Any, ...], dtype[int64]][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: str = 'aspect',
nodeValues=False,
) ndarray[tuple[Any, ...], dtype[floating]][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: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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

Get_SourcePart_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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’

Get_assembly_e(
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

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

Get_columns_e(
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the column indices used to assemble local matrices into the global matrix.

Get_connect_n_e() csr_matrix[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)

Get_dN_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_ddN_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_jacobian_e_pg(
matrixType: MatrixType,
absoluteValues=True,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Returns the jacobians

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

Get_leftDispPart(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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’

Get_list_groupElem(dim=None) list[_GroupElem][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,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

Returns the mesh size of the mesh.

returns meshSize_e if doMean else meshSize_e_s

Get_nPg(matrixType: MatrixType) int[source]#

Returns integration points according to the matrix type.

Get_normals(
nodes: ndarray[tuple[Any, ...], dtype[int64]] | None = None,
displacementMatrix: ndarray[tuple[Any, ...], dtype[floating]] | None = None,
) tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[int64]]][source]#

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

returns normals, nodes.

Get_rows_e(
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the row indices used to assemble local matrices into the global matrix.

Get_weight_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

Returns integration points according to the matrix type.

Get_weightedJacobian_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Returns jacobian_e_pg * weight_pg.

Locates_sol_e(
sol: ndarray[tuple[Any, ...], dtype[floating]],
dof_n: int | None = None,
asFeArray=False,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Locates solution on elements.

Nodes_Circle(
*circles: Circle,
onlyOnEdge=True,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the nodes in the circle(s).

Nodes_Conditions(
func,
) ndarray[tuple[Any, ...], dtype[int64]][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(
*circles: Circle,
direction=[0, 0, 1],
onlyOnEdge=False,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the nodes in the cylinder.

Nodes_Domain(
*domains: Domain,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes in the domain(s).

Nodes_Line(
*lines: Line,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the nodes on the line(s).

Nodes_Point(
*points: Point | ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float],
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes on the point(s).

Nodes_Tags(
*tags: str,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes associated with the tags.

Rotate(
theta: float,
center: tuple = (0, 0, 0),
direction: tuple = (0, 0, 1),
) None[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)

Set_Tag(
nodes: ndarray[tuple[Any, ...], dtype[int64]],
tag: str,
)[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)) None[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)

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

Translates the mesh coordinates.

_Gather() Mesh | None[source]#

Gather the distributed mesh to rank root, reconstructing the global Mesh.

Return type:

Mesh on rank root, None on other ranks.

_Get_realistic_vector_magnitude(coef=0.1) float[source]#

Returns a realistic vector magnitude based on the mesh size.

Parameters:

coef (float, optional) – coef used to scale the average distance between the coordinates and the center, by default 0.1

_ResetMatrix() None[source]#

Resets matrices for each groupElem

copy()[source]#
property Ne: int[source]#

number of elements in the mesh

property Nn: int[source]#

number of nodes in the mesh

property area: float[source]#

total area of the mesh.

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

center of mass / barycenter / inertia center

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

connectivity matrix (Ne, nPe)

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

global nodes coordinates matrix (Nn, 3)

Contains all nodes coordinates

property dict_groupElem: dict[ElemType, _GroupElem][source]#

dictionary containing all the element groups in the mesh

property dim[source]#

mesh dimension

property elemType: ElemType[source]#

mesh element type

property groupElem: _GroupElem[source]#

main group element

property inDim[source]#

dimension in which the mesh lies.

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

property length: float[source]#

total length of the mesh.

property nPe: int[source]#

nodes per element

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

mesh nodes

property orphanNodes: list[int][source]#

nodes not connected to any mesh elements

property verbosity: bool[source]#

the mesh can write in the terminal

property volume: float[source]#

total volume of the mesh.

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

Bases: object

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

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

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

static _Construct_2D_meshes(L=10, h=10, meshSize=3) list[Mesh][source]#

Creates 2D meshes.

static _Construct_3D_meshes(
L=130,
h=13,
b=13,
meshSize=7.5,
useImport3D=False,
) list[Mesh][source]#

Creates 3D meshes.

static _Set_mesh_order(elemType: ElemType) None[source]#

Sets the mesh order

Create_posFile(
coord: ndarray[tuple[Any, ...], dtype[floating]],
values: ndarray[tuple[Any, ...], dtype[floating]],
folder: str,
filename='data',
) str[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

Mesh_2D(
contour: _Geom | Circle | Domain | Points | Contour,
inclusions: list[_Geom | Circle | Domain | Points | Contour] = [],
elemType: ElemType = ElemType.TRI3,
cracks: list[Line | Points | Contour | CircleArc] = [],
refineGeoms: list[Point | Circle | str] = [],
isOrganised=False,
additionalSurfaces: list[_Geom | Circle | Domain | Points | Contour] = [],
additionalLines: list[Line | CircleArc] = [],
additionalPoints: list[Point] = [],
path='',
) Mesh[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.

  • path (str, optional) – path used to save the meshfile, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Beams(
beams: list[_Beam],
elemType=ElemType.SEG2,
additionalPoints: list[Point] = [],
path: str = '',
) Mesh[source]#

Creates a beam mesh.

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

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

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

  • path (str, optional) – path used to save the meshfile, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Extrude(
contour: _Geom | Circle | Domain | Points | Contour,
inclusions: list[_Geom | Circle | Domain | Points | Contour] = [],
extrude: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 0, 1),
layers: list[int] = [],
elemType: ElemType = ElemType.TETRA4,
cracks: list[Line | Points | Contour | CircleArc] = [],
refineGeoms: list[Point | Circle | str] = [],
isOrganised=False,
additionalSurfaces: list[_Geom | Circle | Domain | Points | Contour] = [],
additionalLines: list[Line | CircleArc] = [],
additionalPoints: list[Point] = [],
path='',
) Mesh[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.

  • path (str, optional) – path used to save the meshfile, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Import_mesh(
mesh: str,
setPhysicalGroups=False,
coef=1.0,
) Mesh[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: str,
dim: int,
meshSize=0.0,
elemType: ElemType | None = None,
refineGeoms=[None],
path='',
) Mesh[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

  • path (str, optional) – path used to save the meshfile, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Mesh_Revolve(contour: ~EasyFEA.Geoms._geom._Geom | ~EasyFEA.Geoms._circle.Circle | ~EasyFEA.Geoms._domain.Domain | ~EasyFEA.Geoms._points.Points | ~EasyFEA.Geoms._contour.Contour, inclusions: list[~EasyFEA.Geoms._geom._Geom | ~EasyFEA.Geoms._circle.Circle | ~EasyFEA.Geoms._domain.Domain | ~EasyFEA.Geoms._points.Points | ~EasyFEA.Geoms._contour.Contour] = [], axis: ~EasyFEA.Geoms._line.Line = <EasyFEA.Geoms._line.Line object>, angle=360, layers: list[int] = [30], elemType: ~EasyFEA.FEM._utils.ElemType = ElemType.TETRA4, cracks: list[~EasyFEA.Geoms._line.Line | ~EasyFEA.Geoms._points.Points | ~EasyFEA.Geoms._contour.Contour | ~EasyFEA.Geoms._circle.CircleArc] = [], refineGeoms: list[~EasyFEA.Geoms._utils.Point | ~EasyFEA.Geoms._circle.Circle | str] = [], isOrganised=False, additionalSurfaces: list[~EasyFEA.Geoms._geom._Geom | ~EasyFEA.Geoms._circle.Circle | ~EasyFEA.Geoms._domain.Domain | ~EasyFEA.Geoms._points.Points | ~EasyFEA.Geoms._contour.Contour] = [], additionalLines: list[~EasyFEA.Geoms._line.Line | ~EasyFEA.Geoms._circle.CircleArc] = [], additionalPoints: list[~EasyFEA.Geoms._utils.Point] = [], path='') Mesh[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.

  • path (str, optional) – path used to save the meshfile, by default “” does not save the mesh

Returns:

Created mesh

Return type:

Mesh

Save_Simu(
simu: _Simu,
results: list[str] = [],
details: bool = False,
edgeColor: str = 'black',
plotMesh: bool = True,
showAxes: bool = False,
folder: str = '',
) None[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 “”

Set_meshSize(meshSize: float) None[source]#

Sets the mesh size

_(crack: Contour)[source]#
_Additional_Lines(
dim: int,
lines: list[Line | CircleArc],
) None[source]#

Adds lines to existing dim entities. WARNING: lines must be within the domain.

Parameters:
  • dim (int) – dimension (dim >= 1)

  • lines (list[Union[Line, CircleArc]]) – lines

_Additional_Points(
dim: int,
points: list[Point],
meshSize: float = 0.0,
) None[source]#

Adds points to existing dim entities. WARNING points must be within the domain.

Parameters:
  • dim (int) – dimension (dim >= 1)

  • points (list[Point]) – points

  • meshSize (float) – meshSize

_Additional_Surfaces(
dim: int,
surfaces: list[_Geom | Circle | Domain | Points | Contour],
elemType: ElemType,
isOrganised: bool,
) None[source]#

Adds surfaces to existing dim entities. Tip: if the mesh is not well generated, you can also give the inclusions.

Parameters:
  • dim (int) – dimension (dim >= 2)

  • surfaces (list[Domain | Circle | Points | Contour]) – surfaces

  • elemType (ElemType) – element type used

  • isOrganised (bool) – mesh is organized

_Cracks_SetPhysicalGroups(
cracks: list[Line | Points | Contour | CircleArc],
entities: list[tuple],
) tuple[int | None, int | None, int | None, int | None][source]#

Creates physical groups for cracks embeded in existing gmsh entities.

returns crackLines, crackSurfaces, openPoints, openLines

_Create_Circle(
circle: Circle,
) tuple[int, list[int], list[int]][source]#

Creates a loop with a circle object.

returns loop, lines, points

_Create_Contour(
contour: Contour,
) tuple[int, list[int], list[int], list[int], list[int]][source]#

Creates a loop with a contour object (list of Line, CircleArc, Points).

returns loop, lines, points, openLines, openPoints

_Create_Domain(
domain: Domain,
) tuple[int, list[int], list[int]][source]#

Creates a loop with a domain object.

returns loop, lines, points

_Create_Lines(geom: _Geom, p1, p2) list[int][source]#
_Create_Lines(line: Line, p1, p2)
_Create_Lines(circleArc: CircleArc, p1, p2)
_Create_Lines(points: Points, p1, p2)

Creates the lines in order to construct the contour object (Line, CircleArc, Points).

returns lines

_Extrude(
surfaces: list[int],
extrude: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 0, 1),
elemType: ElemType = ElemType.TETRA4,
layers: list[int] = [],
) list[tuple[int, int]][source]#

Extrudes gmsh surfaces and returns extruded entities.

Parameters:
  • surfaces (list[int]) – gmsh surfaces

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

  • elemType (ElemType, optional) – element type used, by default “TETRA4”

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

_Init_gmsh(factory: str = 'occ') None[source]#

Initializes gmsh.

Links 2 contours and create a volume.

Contours must be connectable, i.e. they must have the same number of points and lines.

Parameters:
  • contour1 (Contour) – the first contour

  • contour2 (Contour) – the second contour

  • elemType (ElemType) – element type used

  • nLayers (int, optional) – number of layers joining the contours, by default 0

  • numElems (list[int], optional) – number of elements for each lines in contour, by default []

Returns:

created entities

Return type:

list[tuple[int, int]]

_Linking_Surfaces(
lines1: list[int],
points1: list[int],
lines2: list[int],
points2: list[int],
linkingLines: list[int],
elemType: ElemType,
nLayers: int = 0,
listPoints: list[list[Point]] = [],
)[source]#

Creates linking surfaces based on lines and linkingLines.

return linkingSurfaces.

_Loop_From_Geom(
geom: _Geom,
) tuple[int, list[int], list[int]][source]#
_Loop_From_Geom(domain: Domain)
_Loop_From_Geom(circle: Circle)
_Loop_From_Geom(points: Points)
_Loop_From_Geom(contour: Contour)

Creates a loop based on the geometric object.

returns loop, lines, points

_Mesh_Generate(
dim: int,
elemType: ElemType,
crackLines: int | None = None,
crackSurfaces: int | None = None,
openPoints: int | None = None,
openLines: int | None = None,
path: str = '',
) None[source]#

Generates the mesh with the available created entities.

Parameters:
  • dim (int) – mesh dimension (e.g 1, 2, 3)

  • elemType (ElemType) – element type

  • crackLines (int, optional) – physical group for crack lines (associated with openPoints), by default None

  • crackSurfaces (int, optional) – physical group for crack surfaces (associated with openLines), by default None

  • openPoints (int, optional) – physical group for open points, by default None

  • openLines (int, optional) – physical group for open lines, by default None

  • path (str, optional) – path used to save the meshfile, by default “” does not save the mesh

_Mesh_Get_Mesh(coef: float = 1.0) Mesh[source]#

Creates the mesh object from the created gmsh mesh.

_Mesh_Refine(
refineGeoms: list[Point | Circle | str],
meshSize: float,
extrude=[0, 0, 1],
) None[source]#

Sets a background mesh

Parameters:
  • refineGeoms (list[Domain|Circle|str]) – Geometric objects to refine de background mesh

  • meshSize (float) – size of elements outside the domain

  • extrude (list[float]) – extrusion vector, by default [0,0,1]

_Revolve(
surfaces: list[int],
axis: Line,
angle: float = 360.0,
elemType: ElemType = ElemType.TETRA4,
layers: list[int] = [30],
) list[tuple[int, int]][source]#

Revolves gmsh surfaces and returns revolved entities.

Parameters:
  • surfaces (list[int]) – gmsh surfaces

  • axis (Line) – rotation axis

  • angle (float, optional) – rotation angle in deg, by default 360.0

  • elemType (ElemType, optional) – element type used

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

_Set_PhysicalGroups(
setPoints: bool = True,
setLines: bool = True,
setSurfaces: bool = True,
setVolumes: bool = True,
) None[source]#

Creates physical groups based on created entities.

_Set_mesh_algorithm(elemType: ElemType) None[source]#

Sets the mesh algorithm.

Mesh.Algorithm

2D mesh algorithm (1: MeshAdapt, 2: Automatic, 3: Initial mesh only, 5: Delaunay, 6: Frontal-Delaunay, 7: BAMG, 8: Frontal-Delaunay for Quads, 9: Packing of Parallelograms, 11: Quasi-structured Quad) Default value: 6

Mesh.Algorithm3D

3D mesh algorithm (1: Delaunay, 3: Initial mesh only, 4: Frontal, 7: MMG3D, 9: R-tree, 10: HXT) Default value: 1

Mesh.RecombinationAlgorithm

Mesh recombination algorithm (0: simple, 1: blossom, 2: simple full-quad, 3: blossom full-quad) Default value: 1

Mesh.SubdivisionAlgorithm

Mesh subdivision algorithm (0: none, 1: all quadrangles, 2: all hexahedra, 3: barycentric) Default value: 0

_Spline_From_Points(
points: Points,
) tuple[int, list[int]][source]#

Creates a gmsh spline from points.

returns spline, points

_Surface_From_Loops(loops: list[int]) int[source]#

Creates a gmsh surface with a gmsh loop (must be a plane surface).

returns surface

_Surfaces(
contour: _Geom | Circle | Domain | Points | Contour,
inclusions: list[_Geom | Circle | Domain | Points | Contour] = [],
elemType: ElemType = ElemType.TRI3,
isOrganised: bool = False,
) tuple[list[int], list[int], list[int]][source]#

Creates gmsh surfaces.

They must be plane surfaces otherwise you must use ‘factory.addSurfaceFilling’ function.

returns surfaces, lines, points

Parameters:
  • contour (Domain | Circle | Points | Contour) – the object that creates the surface area

  • inclusions (list[Domain | Circle | Points | Contour]) –

    hollow or filled objects contained in the contour surface.

    CAUTION: all inclusions must be contained within the contour and must not intersect.

  • elemType (ElemType, optional) – element type used, by default TRI3

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

_Surfaces_Organize(
surfaces: list[int],
elemType: ElemType,
isOrganised: bool = False,
numElems: Collection[int | float] = [],
) None[source]#

Organizes surfaces.

Parameters:
  • surfaces (list[int]) – list of gmsh surfaces

  • elemType (ElemType) – element type

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

  • numElems (_types.Numbers, optional) – number of elements per line, by default []

_Synchronize() None[source]#

Synchronizes geometric entities created on gmsh.

class EasyFEA.FEM._Form(
form: Callable[[...], FeArray | ndarray[tuple[Any, ...], dtype[Any]]],
)[source]#

Bases: ABC

Form class from which BilinearForm and LinearForm are derived.

abstractmethod Assemble(field: Field) csr_matrix[source]#
abstractmethod Integrate_e(field: Field) ndarray[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._GroupElem(
gmshId: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: ABC

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

static _Eval_Functions(
functions: ndarray[tuple[Any, ...], dtype[floating]],
gaussPoints: ndarray[tuple[Any, ...], dtype[floating]],
) ndarray[tuple[Any, ...], dtype[floating]][source]#

Evaluates functions at coordinates.

Use this function to evaluate shape functions.

Parameters:
  • functions (_types.FloatArray) – functions to evaluate, (nPe, dim)

  • gaussPoints (_types.FloatArray) –

    gauss coordinates where functions will be evaluated (nPg, dim).

    Be careful dim must be consistent with function arguments

Returns:

Evaluated functions (nPg, dim, nPe)

Return type:

_types.FloatArray

static _Get_assembly_e(
connect: ndarray[tuple[Any, ...], dtype[int64]],
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

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

Parameters:
  • connect (_types.IntArray) – connectivity matrix (Ne, nPe)

  • dof_n (int) – degree of freedom per node

Get_B_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Get the matrix used to calculate deformations from displacements.

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: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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’

Get_Elements_Nodes(
nodes: ndarray[tuple[Any, ...], dtype[int64]],
exclusively=True,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns elements that exclusively or not use the specified nodes.

Get_Elements_Tag(
tag: str,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns elements associated with the tag.

Get_EulerBernoulli_B_e_pg(
beamStructure: BeamStructure,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Get Euler-Bernoulli beam shape functions derivatives

Get_EulerBernoulli_N_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_EulerBernoulli_N_e_pg_for_beam(
beamStructure: BeamStructure,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Euler-Bernoulli beam shape functions.

Get_EulerBernoulli_N_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[phi_i psi_i … phi_n psi_n]

(nPg, nPe*2)

Get_EulerBernoulli_dN_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_EulerBernoulli_dN_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

(nPg, nPe*2)

Get_EulerBernoulli_ddN_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_EulerBernoulli_ddN_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

(nPg, nPe*2)

Get_F_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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.

Get_GaussCoordinates_e_pg(
matrixType: MatrixType,
elements=array([], dtype=float64),
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

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

Get_Gradient_e_pg(
u: ndarray[tuple[Any, ...], dtype[floating]],
matrixType=MatrixType.rigi,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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

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() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Get_Mapping(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
elements_e: ndarray[tuple[Any, ...], dtype[int64]] | None = None,
needCoordinates=False,
) tuple[ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[floating]] | None][source]#

Locates coordinates within elements.

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_n (The 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_n (The 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: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[N1, … , Nn]

(nPg, 1, nPe)

Get_N_pg_rep(
matrixType: MatrixType,
repeat=1,
) ndarray[tuple[Any, ...], dtype[floating]][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]

Get_Nodes_Circle(
circle: Circle,
onlyOnEdge=False,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes in the circle.

Get_Nodes_Conditions(
func: Callable,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes that meet the specified conditions.

Parameters:

func

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: Circle,
direction=[0, 0, 1],
onlyOnEdge=False,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes in the cylinder.

Get_Nodes_Domain(
domain: Domain,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes in the domain.

Get_Nodes_Line(
line: Line,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes on the line.

Get_Nodes_Point(
point: Point | ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float],
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns nodes on the point.

Get_Nodes_Tag(
tag: str,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns node associated with the tag.

Get_ReactionPart_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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

Get_SourcePart_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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’

Get_assembly_e(
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

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

Parameters:

dof_n (int) – degree of freedom per node

Get_columns_e(
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the column indices used to assemble local matrices into the global matrix.

Get_connect_n_e() csr_matrix[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)

Get_dN_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_dN_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ … Nn,ξ

Ni,η … Nn,η

Ni,ζ … Nn,ζ

(nPg, dim, nPe)

Get_ddN_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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)

Get_ddN_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[Ni,ξ2 … Nn,ξ2

Ni,η2 … Nn,η2

Ni,ζ2 … Nn,ζ2]

(nPg, dim, nPe)

Get_dddN_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[Ni,ξ3 … Nn,ξ3

Ni,η3 … Nn,η3

Ni,ζ3 … Nn,ζ3]

(nPg, dim, nPe)

Get_ddddN_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[Ni,ξ4 … Nn,ξ4

Ni,η4 … Nn,η4

Ni,ζ4 … Nn,ζ4]

(pg, dim, nPe)

Get_gauss(
matrixType: MatrixType,
) Gauss[source]#

Returns integration points according to the matrix type.

Get_invF_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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

Get_jacobian_e_pg(
matrixType: MatrixType,
absoluteValues=True,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Returns the jacobians.

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

Get_leftDispPart(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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’

Get_normals_e_pg(
matrixType: MatrixType,
displacementMatrix: ndarray[tuple[Any, ...], dtype[floating]] | None = None,
) FeArray[source]#

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

Get_pointsInElem(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
elem: int,
) ndarray[tuple[Any, ...], dtype[int64]][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_rows_e(
dof_n: int,
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Returns the row indices used to assemble local matrices into the global matrix.

Get_weight_pg(
matrixType: MatrixType,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

Returns integration point weights according to the matrix type.

Get_weightedJacobian_e_pg(
matrixType: MatrixType,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Returns the jacobian_e_pg * weight_pg.

Integrate_e(
func=<function _GroupElem.<lambda>>,
matrixType=MatrixType.mass,
) ndarray[tuple[Any, ...], dtype[floating]][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: ndarray[tuple[Any, ...], dtype[floating]],
dof_n: int | None = None,
asFeArray=False,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Locates sol on elements

Set_Tag(
nodes: ndarray[tuple[Any, ...], dtype[int64]],
tag: str,
)[source]#

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

_EulerBernoulli_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

_EulerBernoulli_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

(nPe*2, 1)

_EulerBernoulli_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

(nPe*2, 2)

_Get_Mapping(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
elements_e: ndarray[tuple[Any, ...], dtype[int64]],
needCoordinates=False,
) tuple[ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[floating]] | None][source]#

Locates coordinates within elements.

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_n (The 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_n (The coordinates of the identified nodes, expressed in the reference element’s (ξ, η, ζ) coordinate system.) – This is applicable only if needCoordinates is set to True.

_Get_coord_Near(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
coordElem: ndarray[tuple[Any, ...], dtype[floating]],
dims: ndarray[tuple[Any, ...], dtype[floating]],
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Get indexes in coordinates_n that are within the coordElem’s bounds.

Parameters:
  • coordinates_n (_types.FloatArray) – coordinates to check

  • coordElem (_types.FloatArray) – element’s bounds

  • dims (_types.FloatArray) – (nX, nY, nZ) = np.max(coordinates_n, 0) - np.min(coordinates_n, 0) + 1

Returns:

indexes in element’s bounds.

Return type:

_types.IntArray

_Get_nearby_elements(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Get nearby elements.

Parameters:

coordinates_n (_types.FloatArray) – coordinates (n, 3) array

Returns:

nearby elements

Return type:

_types.IntArray

_Get_nearby_nodes(
coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
) ndarray[tuple[Any, ...], dtype[int64]][source]#

Get nearby nodes.

Parameters:

coordinates_n (_types.FloatArray) – coordinates (n, 3) array

Returns:

nearby nodes

Return type:

_types.IntArray

_Get_partitioned_data() tuple[ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[int64]], ndarray[tuple[Any, ...], dtype[int64]]][source]#

Returns the partitioned data used in mpi.

(rank, elements, ghostElements, nodes, ghostNodes)

_Get_sysCoord_e(
displacementMatrix: ndarray[tuple[Any, ...], dtype[Any]] | None = None,
)[source]#

Get the basis transformation matrix (Ne, 3, 3).

[ix, jx, kx

iy, jy, ky

iz, jz, kz]

This matrix can be used to project points with (x, y, z) coordinates into the element’s (i, j, k) coordinate system.

coordo_e * sysCoord_e -> coordinates in element’s (i, j, k) coordinate system.

_InitMatrix() None[source]#

Initializes matrix dictionaries for finite element construction

_Init_Functions(
order: int,
) ndarray[tuple[Any, ...], dtype[floating]][source]#

Initializes functions to be evaluated at Gauss points.

abstractmethod _N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_Set_partitioned_data(
elements: ndarray[tuple[Any, ...], dtype[int64]],
nodes: ndarray[tuple[Any, ...], dtype[int64]],
rank: int = 0,
ghostElements: ndarray[tuple[Any, ...], dtype[int64]] = [],
) None[source]#

Sets the partitioned data used in mpi.

Parameters:
  • elements (_types.IntArray) – the positions of elements in the global mesh.

  • nodes (_types.IntArray) – the (non-ghost) nodes.

  • rank (int, optional) – mpi rank, by default 0

  • ghostElements (_types.IntArray, optional) – the positions of ghost elements in the global mesh, by default [].

  • Remark

  • ------

  • array. (Ghost nodes will be computed using the given (non-ghost) nodes)

abstractmethod _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

abstractmethod _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

abstractmethod _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

abstractmethod _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

property Ne: int[source]#

number of elements

property Nedge: int[source]#

number of edge nodes per element

property Nface: int[source]#

number of face nodes per element

property Nn: int[source]#

number of nodes used by the element group

property Nvertex: int[source]#

number of vertex nodes per element

property Nvolume: int[source]#

number of volume nodes per element

property _dict_elements_tags: dict[str, ndarray[tuple[Any, ...], dtype[int64]]][source]#

dictionary associating tags with elements.

property _dict_nodes_tags: dict[str, ndarray[tuple[Any, ...], dtype[int64]]][source]#

Dictionary associating tags with nodes.

property _rankPartition: ndarray[tuple[Any, ...], dtype[int64]] | None[source]#

Array of shape (Ne,) where rankPartition[e] is the MPI rank that owned element e. Set on rank 0 only after Gather(). None on all other ranks.

property area: float[source]#

area covered by elements

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

area covered by each element

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

center of mass / barycenter / inertia center

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

connectivity matrix (Ne, nPe)

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

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

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

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

property dim: int[source]#

element dimension

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

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

property elemType: ElemType[source]#

element type

property elementTags: list[str][source]#

returns element tags.

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

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

property gmshId: int[source]#

gmsh Id

property inDim: int[source]#

dimension in which the elements are located

property length: float[source]#

length covered by elements

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

length covered by each element

property nPe: int[source]#

nodes per element

property nodeTags: list[str][source]#

Returns node tags.

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

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

property order: int[source]#

element order

abstract property origin: list[int][source]#

reference element origin coordinates

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

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

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

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[source]#

element topology

abstract property triangles: list[int][source]#

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

property volume: float[source]#

volume covered by elements

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

volume covered by each element

EasyFEA.FEM.Calc_projector(
oldMesh: Mesh,
newMesh: Mesh,
) csr_matrix[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: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Computes det(mat)

EasyFEA.FEM.Inv(
mat: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
)[source]#

Computes inv(mat)

EasyFEA.FEM.Mesh_Optim(
DoMesh: Callable[[str], Mesh],
folder: str,
criteria: str = 'aspect',
quality=0.8,
ratio: float = 0.7,
iterMax=20,
coef: float = 0.5,
) tuple[Mesh, float][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]

EasyFEA.FEM.Norm(
array: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
**kwargs,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

np.linalg.norm() wrapper.

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

EasyFEA.FEM.Sym_Grad(
u: Field,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
EasyFEA.FEM.TensorProd(
A: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
B: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
symmetric=False,
ndim: int | None = None,
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][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: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Computes trace(mat)

EasyFEA.FEM.Transpose(
mat: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
) FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#

Computes transpose(mat)

class EasyFEA.FEM.Elems.HEXA20(
gmshId: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_EulerBernoulli_N() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

_EulerBernoulli_dN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 1)

_EulerBernoulli_ddN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 2)

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_EulerBernoulli_N() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

_EulerBernoulli_dN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 1)

_EulerBernoulli_ddN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 2)

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_EulerBernoulli_N() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

_EulerBernoulli_dN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 1)

_EulerBernoulli_ddN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 2)

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_EulerBernoulli_N() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

[phi_i psi_i … phi_n psi_n]

(nPe*2, 1)

_EulerBernoulli_dN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 1)

_EulerBernoulli_ddN() ndarray[tuple[Any, ...], dtype[Any]][source]#

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

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

(nPe*2, 2)

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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: int,
connect: ndarray[tuple[Any, ...], dtype[int64]],
coordGlob: ndarray[tuple[Any, ...], dtype[floating]],
)[source]#

Bases: _GroupElem

Get_Local_Coords()[source]#

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

_N() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

N1

Nn

(nPe, 1)

_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

Ni,ξ Ni,η Ni,ζ

Nn,ξ Nn,η Nn,ζ

(nPe, dim)

_ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

_ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#

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

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

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

(nPe, dim)

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

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

property origin: list[int][source]#

reference element origin coordinates

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

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][source]#

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