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
Meshes examples
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:
Mesh_2D(): Generates a 2D mesh.Mesh_Extrude(): Creates a mesh by extruding a 2D shape.Mesh_Revolve(): Generates a mesh by revolving a 2D shape around an axis.Mesh_Import_part(): Imports a CAD part (e.g., .stp) to create a mesh.Mesh_Import_mesh(): Imports an existing Gmsh mesh. EasyFEA is also linked to meshio and can be used through the following functions:Medit_to_EasyFEA(): Imports a Medit mesh.Gmsh_to_EasyFEA(): Imports a Gmsh mesh.PyVista_to_EasyFEA(): Imports a PyVista mesh (UnstructuredGrid or MultiBlock).Ensight_to_EasyFEA(): Imports an EnSight mesh.
See also
Operators#
The EasyFEA.FEM.Operators module provides the element-level operators that integrate a form over the Gauss points and produce the element matrices/vectors assembled into the global system (see Understand the solve pipeline). Below, \(c\) is a scalar/field coefficient, \(\Brm\) the strain-displacement operator and \(\Nrm\) the shape functions.
Bilinear — element matrices#
UV()— \(\int_\Omega c \, u \, v \, \dO\) (mass).GradUGradV()— \(\int_\Omega c \, \grad u \cdot \grad v \, \dO\) (diffusion / Laplacian).GradU_A_GradV()— \(\int_\Omega c \, \grad u \cdot \Abf \cdot \grad v \, \dO\) (anisotropic diffusion).LinearizedElasticity()— \(\int_\Omega \Eps(u) : \Cbf : \Eps(v) \, \dO\) (small-strain stiffness).MassAlongNormal()— \(\int_\Gamma c \, (u \cdot \nbf)(v \cdot \nbf) \, \dS\) (surface mass projected on the normal; Robin terms).BeamBending()— \(\int_e \Brm^\top \Dbf_{\text{bend}} \, \Brm \, dx\) (axial + bending + torsion).BeamShear()— \(\int_e \Brm^\top \Dbf_{\text{shear}} \, \Brm \, dx\) (transverse shear, selective reduced integration).BeamStiffness()— \(\Krm_e = \text{BeamBending} + \text{BeamShear}\).BeamMass()— \(\int_e c \, \Nrm^\top \Mbf \, \Nrm \, dx\) (consistent beam mass).
Linear — element load vectors#
V()— \(\int_\Omega \fbf \cdot v \, \dO\) (body / surface load).
NonLinear — tangent + residual#
SecondPiolaKirchhoffStressTensor()— residual \(\Frm_e = \int \Brm^\top \boldsymbol{\Sigma} \, \dO\) and consistent tangent \(\Krm_e = \int \Brm^\top \dpartial{\boldsymbol{\Sigma}}{\Erm} \Brm \, \dO + \int \grad^\top \boldsymbol{\Sigma} \, \grad \, \dO\) (material + geometric), with \(\boldsymbol{\Sigma}\) the second Piola-Kirchhoff stress.KelvinVoigtDamping()— damping \(\Crm_e = c \, \eta \int \Brm^\top \Brm \, \dO\) and the configuration tangent \(\Krm_e^{\text{geo}} = \dpartial{(\Crm \vrm)}{\urm}\) (large-strain Kelvin–Voigt).FollowingPressure()— follower-pressure load \(\int_\Gamma p \, v \cdot \nbf(u) \, \dS\) with a deformation-dependent normal, plus its (non-symmetric) tangent.
FEM API#
- class EasyFEA.FEM.BiLinearForm(
- form: Callable[[...], FeArray | ndarray[tuple[Any, ...], dtype[Any]]],
Bases:
_FormBilinear form.
- 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,
Bases:
objectBoundary condition.
- static Get_dofs(
- problemType: str,
- list_Bc_Condition: list[BoundaryCondition],
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],
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],
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],
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
- class EasyFEA.FEM.ElemType(value)[source]#
Bases:
str,EnumImplemented Lagrange isoparametric 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'#
- 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,
- static broadcast(
- value,
- Ne: int,
- nPg: int,
- tensor_ndim: int = 0,
Broadcast a scalar or array coefficient to a shape compatible with multiplication against an
(Ne, nPg, ...)FeArray.Returns a stride-tricked read-only view for non-scalar inputs (no data duplication). Callers must not mutate the result in place; the expected use is consumption inside expressions such as
coef * wJ_e_pg * dN_e_pg.T @ dN_e_pg, which create new arrays.tensor_ndimdeclares how many trailing axes ofvalueare tensor dims (e.g.2for a Hooke tensor(..., nstrain, nstrain)). With it set, the leading axes are checked against()/(Ne,)/(Ne, nPg)strictly — this disambiguates shapes like(Ne, n, n)from(Ne, nPg, n)whennPg == n(e.g. TRI6 in 2D, wherenPg == nstrain == 3).Accepted shapes (with default
tensor_ndim=0)#scalar (int / float / numpy scalar) → returned as
float.(Ne, nPg, ...)ndarray / FeArray → wrapped as FeArray.1-D
(Ne,)or(nPg,)→ tiled to(Ne, nPg).Any other shape broadcastable to
(Ne, nPg, ...)→ tiled with leading(Ne, nPg)dims.
- _asfearrays(
- *,
- broadcastFeArrays=False,
- dot(other, /, out=None)[source]#
Refer to
numpy.dot()for full documentation.See also
numpy.dotequivalent function
- integrate() ndarray[source]#
Integrate over the Gauss-point axis (axis 1). Returns
(Ne, ...)ndarray.
- 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
- class EasyFEA.FEM.Field(
- groupElem: _GroupElem,
- dof_n: int,
- matrixType: MatrixType = MatrixType.mass,
Bases:
objectField class.
- Evaluate_e( ) 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.
- 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:
- 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,
Bases:
objectGauss quadrature.
- static Gauss_factory(
- elemType: ElemType,
- matrixType: MatrixType,
Returns the integration points and weights based on the element and matrix type.
- static _Gauss_factory_nPg(
- elemType: ElemType,
- nPg,
Returns the integration points and weights based on the element type and the number of Gauss points (nPg).
- static _Hexahedron(
- nPg: int,
available [8, 27]
order = [3, 5]
- static _Prism(
- nPg: int,
available [6, 8, 21]
order X = [3, 3, 5]
order Y & Z = [2, 3, 5]
- static _Quadrangle(
- nPg: int,
available [4, 9]
order = [1, 2]
- static _Tetrahedron(
- nPg: int,
available [1, 4, 5, 15]
order = [1, 2, 3, 5]
- class EasyFEA.FEM.GroupElemFactory[source]#
Bases:
object- static Create(
- elemType: ElemType,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Creates an element group
- Parameters:
elemType (ElemType) – element type
connect (_types.IntArray) – connection matrix storing nodes for each element (Ne, nPe)
coordinates (_types.FloatArray) – coordinate matrix (contains all mesh coordinates)
- Returns:
the element group
- Return type:
GroupeElem
- static Get_ElemInFos(
- gmshId: int,
return elemType, nPe, dim, order, Nvertex, Nedge, Nface, Nvolume
associated with the gmsh id.
- static _Create(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Creates an element group.
- Parameters:
gmshId (int) – id gmsh
connect (_types.IntArray) – connection matrix storing nodes for each element (Ne, nPe)
coordinates (_types.FloatArray) – coordinate matrix (contains all mesh coordinates)
- Returns:
the element group
- Return type:
GroupeElem
- static _Get_2d_element_types(
- elemType: ElemType,
Returns 2d element types associated with the elementType.
- 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:
- 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:
- GROUP_CLASS_MAP: dict[ElemType, _GroupElem] = {ElemType.HEXA20: <class 'EasyFEA.FEM.Elems._hexa.HEXA20'>, ElemType.HEXA27: <class 'EasyFEA.FEM.Elems._hexa.HEXA27'>, ElemType.HEXA8: <class 'EasyFEA.FEM.Elems._hexa.HEXA8'>, ElemType.POINT: <class 'EasyFEA.FEM.Elems._point.POINT'>, ElemType.PRISM15: <class 'EasyFEA.FEM.Elems._prism.PRISM15'>, ElemType.PRISM18: <class 'EasyFEA.FEM.Elems._prism.PRISM18'>, ElemType.PRISM6: <class 'EasyFEA.FEM.Elems._prism.PRISM6'>, ElemType.QUAD4: <class 'EasyFEA.FEM.Elems._quad.QUAD4'>, ElemType.QUAD8: <class 'EasyFEA.FEM.Elems._quad.QUAD8'>, ElemType.QUAD9: <class 'EasyFEA.FEM.Elems._quad.QUAD9'>, ElemType.SEG2: <class 'EasyFEA.FEM.Elems._seg.SEG2'>, ElemType.SEG3: <class 'EasyFEA.FEM.Elems._seg.SEG3'>, ElemType.SEG4: <class 'EasyFEA.FEM.Elems._seg.SEG4'>, ElemType.SEG5: <class 'EasyFEA.FEM.Elems._seg.SEG5'>, ElemType.TETRA10: <class 'EasyFEA.FEM.Elems._tetra.TETRA10'>, ElemType.TETRA4: <class 'EasyFEA.FEM.Elems._tetra.TETRA4'>, ElemType.TRI10: <class 'EasyFEA.FEM.Elems._tri.TRI10'>, ElemType.TRI15: <class 'EasyFEA.FEM.Elems._tri.TRI15'>, ElemType.TRI3: <class 'EasyFEA.FEM.Elems._tri.TRI3'>, ElemType.TRI6: <class 'EasyFEA.FEM.Elems._tri.TRI6'>}#
- 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 = '',
Bases:
BoundaryConditionLagrange boundary condition
- class EasyFEA.FEM.LinearForm(
- form: Callable[[...], FeArray | ndarray[tuple[Any, ...], dtype[Any]]],
Bases:
_FormLinear form.
- class EasyFEA.FEM.MatrixType(value)[source]#
Bases:
str,EnumOrder 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
- beam_shear = 'beam_shear'#
Reduced integration order for the Timoshenko shear term (n-1 Gauss points for SEGn). Selective reduced integration is the standard cure for shear locking in Lagrange Timoshenko beams (Hughes 1987; Bathe 2006).
- 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:
ObservableMesh class that contains several _GroupElem instances.
- static Merge(
- list_mesh: list[Mesh],
- constructUniqueElements: bool = True,
- mergePoints: bool = True,
- mergePointsTol: float = 1e-12,
- return_mapping: bool = False,
Merges EasyFEA meshes into a single mesh.
This is a static method: call it as
Mesh.Merge([mesh1, mesh2, ...]).- Parameters:
list_mesh (list[Mesh]) – Meshes to merge.
constructUniqueElements (bool, optional) – Remove duplicate elements after merging, by default True.
mergePoints (bool, optional) – Merge coincident nodes across meshes, by default True. Set to False to skip the KDTree search and simply concatenate coordinate arrays (faster when meshes are known to be disjoint).
mergePointsTol (float, optional) – Absolute distance below which two points are considered identical and merged, by default 1e-12. Ignored when
mergePoints=False.return_mapping (bool, optional) – If True, also return the node mapping as a list of integer arrays, one per input mesh.
mapping[i][j]is the index of nodejfromlist_mesh[i]in the merged mesh, by default False.
- Returns:
Mesh – The merged mesh.
mapping (list[np.ndarray], only when
return_mapping=True) –mapping[i]maps each node index oflist_mesh[i]to its index in the merged mesh. Useful to transfer node tags:merged, mapping = Mesh.Merge([m1, m2], return_mapping=True) for tag in m1.groupElem.nodeTags: nodes = mapping[0][m1.groupElem.Get_Nodes_Tag(tag)] merged.groupElem.Set_Tag(nodes, tag)
- 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,
Returns elements that exclusively or not use the specified nodes.
- Elements_Tags(
- *tags: str,
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,
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_New_meshSize_n(
- error_e: ndarray[tuple[Any, ...], dtype[floating]],
- coef: float = 0.5,
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]],
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,
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,
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_assembly_e(
- dof_n: int,
Returns assembly matrix for specified dof_n (Ne, nPe*dof_n)
- Get_columns_e(
- dof_n: int,
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_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,
Returns the mesh size of the mesh.
returns meshSize_e if doMean else meshSize_e_s
- Get_normals(
- nodes: ndarray[tuple[Any, ...], dtype[int64]] | None = None,
- displacementMatrix: ndarray[tuple[Any, ...], dtype[floating]] | None = None,
Returns normal vectors and nodes belonging to the edge of the mesh.
returns normals, nodes.
- Get_rows_e(
- dof_n: int,
Returns the row indices used to assemble local matrices into the global matrix.
- Locates_sol_e(
- sol: ndarray[tuple[Any, ...], dtype[floating]],
- dof_n: int | None = None,
- asFeArray=False,
Locates solution on elements.
- Nodes_Circle(
- *circles: Circle,
- onlyOnEdge=True,
Returns the nodes in the circle(s).
- Nodes_Conditions(
- func,
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,
Returns the nodes in the cylinder.
- Nodes_Domain(
- *domains: Domain,
Returns nodes in the domain(s).
- Nodes_Line(
- *lines: Line,
Returns the nodes on the line(s).
- Nodes_Point(
- *points: Point | ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float],
Returns nodes on the point(s).
- Nodes_Tags(
- *tags: str,
Returns nodes associated with the tags.
- Rotate(
- theta: float,
- center: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 0, 0),
- direction: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 0, 1),
Rotates the mesh coordinates around an axis defined by a center and a direction.
- Parameters:
theta (float) – rotation angle [deg]
center (_types.Coords, optional) – point on the rotation axis, by default (0,0,0)
direction (_types.Coords, optional) – unit vector defining the rotation axis, by default (0,0,1) (z-axis)
- Save(folder: str, filename: str = 'mesh', gather=False) str[source]#
Saves the mesh in the folder. Saves the mesh as ‘filename.pickle’.
- Set_Tag(
- nodes: ndarray[tuple[Any, ...], dtype[int64]],
- tag: str,
Set a tag on the nodes and elements belonging to each group of elements in the mesh.
- Symmetry(
- point: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 0, 0),
- n: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (1, 0, 0),
Reflects the mesh coordinates through a plane defined by a point and a normal vector.
- Parameters:
point (_types.Coords, optional) – a point on the reflection plane, by default (0,0,0)
n (_types.Coords, optional) – normal vector to the reflection plane, by default (1,0,0) (yz-plane)
- Translate(dx: float = 0.0, dy: float = 0.0, dz: float = 0.0) None[source]#
Translates the mesh coordinates in 3D space.
- Parameters:
dx (float, optional) – translation along the x-axis, by default 0.0
dy (float, optional) – translation along the y-axis, by default 0.0
dz (float, optional) – translation along the z-axis, by default 0.0
- _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
- property center: ndarray[tuple[Any, ...], dtype[floating]][source]#
center of mass / barycenter / inertia center
- property coord: ndarray[tuple[Any, ...], dtype[floating]][source]#
global nodes coordinates matrix (Ncoords, 3)
Contains all nodes coordinates
- property dict_groupElem: dict[ElemType, _GroupElem][source]#
dictionary containing all the element groups in the mesh
- property groupElem: _GroupElem[source]#
main group element
- class EasyFEA.FEM.Mesher(openGmsh: bool = False, gmshVerbosity: bool = False, verbosity: bool = False)[source]#
Bases:
objectMesher class used to construct and generate the mesh via gmsh.
- static Get_Entities(
- points=[],
- lines=[],
- surfaces=[],
- volumes=[],
Get entities from from points, lines, surfaces and volumes tags
- static _Construct_3D_meshes(
- L=130,
- h=13,
- b=13,
- meshSize=7.5,
- useImport3D=False,
Creates 3D meshes.
- Create_posFile(
- coord: ndarray[tuple[Any, ...], dtype[floating]],
- values: ndarray[tuple[Any, ...], dtype[floating]],
- folder: str,
- filename='data',
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='',
Creates a 2D mesh from a contour and inclusions that must form a closed plane surface.
- Parameters:
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_Beams( ) Mesh[source]#
Creates a beam mesh.
- Parameters:
beams (list[_Beam]) – list of Beams
elemType (ElemType, optional) – element type, by default “SEG3” [“SEG2”, “SEG3”, “SEG4”, “SEG5”]
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_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='',
Creates a 3D mesh by extruding a surface constructed from a contour and inclusions.
- Parameters:
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_Import_mesh(
- mesh: str,
- setPhysicalGroupsFromEntities=True,
- meshOrder: int = None,
- coef=1.0,
Creates the mesh from an .msh file.
- Parameters:
mesh (str) – .msh file
setPhysicalGroupsFromEntities (bool, optional) – reconstruct physical groups from gmsh entities, by default True. Each entity becomes a physical group named
f"{P|L|S|V}{tag}"wheretagis the entity’s own gmsh tag (so a volume entity with tag 3 becomes"V3"). Existing physical-group names are preserved.meshOrder (int, optional) – mesh order to use in gmsh.model.mesh.setOrder(meshOrder)`, by default None
coef (float, optional) – coef applied to the node coordinates, by default 1.0
- Returns:
Created mesh
- Return type:
- Mesh_Import_part(
- file: str,
- dim: int,
- meshSize=0.0,
- elemType: ElemType | None = None,
- refineGeoms=[None],
- path='',
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_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:
- Save_Simu(
- simu: _Simu,
- results: list[str] = [],
- details: bool = False,
- edgeColor: str = 'black',
- plotMesh: bool = True,
- showAxes: bool = False,
- folder: str = '',
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 “”
- _Additional_Lines( ) None[source]#
Adds lines to existing dim entities. WARNING: lines must be within the domain.
- _Additional_Points(
- dim: int,
- points: list[Point],
- meshSize: float = 0.0,
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( ) None[source]#
Adds surfaces to existing dim entities. Tip: if the mesh is not well generated, you can also give the inclusions.
- _Cracks_SetPhysicalGroups( ) 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,
Creates a loop with a circle object.
returns loop, lines, points
- _Create_Contour(
- contour: Contour,
Creates a loop with a contour object (list of Line, CircleArc, Points).
returns loop, lines, points, openLines, openPoints
- _Create_Domain(
- domain: Domain,
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
- _Cyclic_Boundary_Order(line_tags: list[int]) list[int][source]#
Returns line_tags in cyclic order around the boundary.
gmsh.model.getBoundaryreturns curves sorted by tag, but transfinite opposite-side pairing needs traversal order. Walk the closed loop by following shared endpoints from one curve to the next.Falls back to the input order if the boundary is not a single closed loop (e.g. multiple holes).
- _Extrude(
- surfaces: list[int],
- extrude: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 0, 1),
- elemType: ElemType = ElemType.TETRA4,
- layers: list[int] = [],
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 []
- _Link_Contours(
- contour1: Contour,
- contour2: Contour,
- elemType: ElemType,
- nLayers: int = 0,
- numElems: list[int] = [],
Links 2 contours and create a volume.
Contours must be connectable, i.e. they must have the same number of points and lines.
- Parameters:
- 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]] = [],
Creates linking surfaces based on lines and linkingLines.
return linkingSurfaces.
- _Loop_From_Geom(
- geom: _Geom,
- _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 = '',
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
- _Organise_Surfaces(
- elemType: ElemType,
- isOrganised: bool,
- meshSize: float,
Applies transfinite settings to every final 2D entity.
For each surface, walks the boundary in cyclic order, computes per-line node counts from the actual line length queried via
gmsh.model.occ.getMassand the givenmeshSize, equalizes opposite-side pairs (transfinite quad requires this), and appliessetTransfiniteCurve/setTransfiniteSurface/setRecombine.meshSize == 0means “no user sizing”, so transfinite is skipped and gmsh falls back to its default mesher (matches the legacy behavior when transfinite curve settings were lost across boolean ops).
- _Revolve(
- surfaces: list[int],
- axis: Line,
- angle: float = 360.0,
- elemType: ElemType = ElemType.TETRA4,
- layers: list[int] = [30],
Revolves gmsh surfaces and returns revolved entities.
- _Set_PhysicalGroups(
- setPoints: bool = True,
- setLines: bool = True,
- setSurfaces: bool = True,
- setVolumes: bool = True,
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,
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] = [],
Creates gmsh surfaces (pure geometry — no transfinite settings).
They must be plane surfaces otherwise you must use ‘factory.addSurfaceFilling’ function.
returns surfaces, lines, points.
- _Surfaces_Organize(
- surfaces: list[int],
- elemType: ElemType,
- isOrganised: bool = False,
- numElems: Collection[int | float] = [],
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 []
- class EasyFEA.FEM._Form(
- form: Callable[[...], FeArray | ndarray[tuple[Any, ...], dtype[Any]]],
Bases:
ABCForm class from which BilinearForm and LinearForm are derived.
- class EasyFEA.FEM._GroupElem(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
ABCThe _GroupElem base class, from which all element types inherit.
- static _Eval_Functions(
- functions: ndarray[tuple[Any, ...], dtype[floating]],
- gaussPoints: ndarray[tuple[Any, ...], dtype[floating]],
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,
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,
Strain-displacement operator
ε = B · uin Kelvin–Mandel notation.With
c = 1/√2(the Kelvin–Mandel scaling on shear terms) and DOF columns ordered(uₓ₁, u_y₁, uₓ₂, u_y₂, …, uₓₙ, u_yₙ)in 2D (extend withu_zin 3D):2D — 3 strain components
(ε_xx, ε_yy, √2·ε_xy):` [ N₁,x 0 N₂,x 0 ⋯ Nₙ,x 0 ] [ 0 N₁,y 0 N₂,y ⋯ 0 Nₙ,y ] [ c·N₁,y c·N₁,x c·N₂,y c·N₂,x ⋯ c·Nₙ,y c·Nₙ,x ] `3D — 6 strain components
(ε_xx, ε_yy, ε_zz, √2·ε_yz, √2·ε_xz, √2·ε_xy):` [ N₁,x 0 0 ⋯ Nₙ,x 0 0 ] [ 0 N₁,y 0 ⋯ 0 Nₙ,y 0 ] [ 0 0 N₁,z ⋯ 0 0 Nₙ,z ] [ 0 c·N₁,z c·N₁,y ⋯ 0 c·Nₙ,z c·Nₙ,y] [ c·N₁,z 0 c·N₁,x ⋯ c·Nₙ,z 0 c·Nₙ,x] [ c·N₁,y c·N₁,x 0 ⋯ c·Nₙ,y c·Nₙ,x 0 ] `The
cfactor makesBorthonormal soεᵀ·σequals the inner productε:σdirectly — no factor-of-2 correction like Voigt needs.Shape:
(Ne, nPg, 3, nPe·dim)in 2D,(Ne, nPg, 6, nPe·dim)in 3D.
- Get_DiffusePart_e_pg(
- matrixType: MatrixType,
Building block of the scalar diffusion term — precomputed
wJ · dNᵀ.` diffusePart_e_pg = wJ · dNᵀ D1_e = (diffusePart_e_pg · A · dN).sum(over Gauss points) D2_e = (diffusePart_e_pg · dN).sum(over Gauss points) `Shape:
(Ne, nPg, nPe, dim).
- Get_Elements_Nodes(
- nodes: ndarray[tuple[Any, ...], dtype[int64]],
- exclusively=True,
Returns elements that exclusively or not use the specified nodes.
- Get_Elements_Tag(
- tag: str,
Returns elements associated with the tag.
- Get_F_e_pg(
- matrixType: MatrixType,
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),
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,
Displacement gradient
∇uat the Gauss points, padded to3×3.- Parameters:
u – Flattened nodal field
[ux₁, uy₁, uz₁, …, uxₙ, uyₙ, uzₙ]of sizeNn · dim.dimis inferred fromu.size.matrixType – Integration scheme. Defaults to
MatrixType.rigi.
- Returns:
Gradient field. Layout depends on
dim(unused components are zeroed so the result is always3×3):dim == 1:` [ ux,x 0 0 ] [ 0 0 0 ] [ 0 0 0 ] `dim == 2:` [ ux,x ux,y 0 ] [ uy,x uy,y 0 ] [ 0 0 0 ] `dim == 3:` [ ux,x ux,y ux,z ] [ uy,x uy,y uy,z ] [ uz,x uz,y uz,z ] `Shape:
(Ne, nPg, 3, 3).- Return type:
- 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,
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,
Shape functions evaluated at Gauss points in (ξ, η, ζ).
` [ N₁ ⋯ Nₙ ] `Shape:(nPg, 1, nPe).
- Get_N_pg_rep(
- matrixType: MatrixType,
- repeat=1,
Shape functions arranged in a block-diagonal matrix of width
repeat.repeat=1— scalar fields (thermal, phase-field):` [ N₁ ⋯ Nₙ ] `repeat=2— 2D vector fields (elastic):` [ N₁ 0 ⋯ Nₙ 0 ] [ 0 N₁ ⋯ 0 Nₙ ] `repeat=3— 3D vector fields: analogous3 × 3·nPeblock-diagonal pattern.Shape:
(nPg, repeat, repeat·nPe).
- Get_Nodes_Circle(
- circle: Circle,
- onlyOnEdge=False,
Returns nodes in the circle.
- Get_Nodes_Conditions(
- func: Callable,
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,
Returns nodes in the cylinder.
- Get_Nodes_Domain(
- domain: Domain,
Returns nodes in the domain.
- Get_Nodes_Line(
- line: Line,
Returns nodes on the line.
- Get_Nodes_Point(
- point: Point | ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float],
Returns nodes on the point.
- Get_Nodes_Tag(
- tag: str,
Returns node associated with the tag.
- Get_ReactionPart_e_pg(
- matrixType: MatrixType,
- dof_n: int = 1,
Building block of the reaction/mass term — precomputed
wJ · Nᵀ · N.dof_n=1(default) — scalar fields (thermal, phase-field).dof_n=dim— vector fields (elastic dynamics); uses the block-diagonal shape functions fromGet_N_pg_rep().` reactionPart_e_pg = wJ · Nᵀ · N M_e = (coef · reactionPart_e_pg).sum(over Gauss points) `Shape:
(Ne, nPg, nPe·dof_n, nPe·dof_n).
- Get_SourcePart_e_pg(
- matrixType: MatrixType,
- dof_n: int = 1,
Building block of the source/body-force term — precomputed
wJ · Nᵀ.dof_n=1(default) — scalar fields (thermal heat source, phase-field).dof_n=dim— vector body forces (elastic, hyperelastic); uses the block-diagonal shape functions fromGet_N_pg_rep().` sourcePart_e_pg = wJ · Nᵀ F_e = (f · sourcePart_e_pg).sum(over Gauss points) `Same pattern as
Get_ReactionPart_e_pg()but for the linear form∫ f · v dΩ. Eachdof_nvalue is cached independently.Shape:
(Ne, nPg, nPe·dof_n, dof_n).
- Get_assembly_e(
- dof_n: int,
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,
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,
First derivatives of the shape functions in physical (x, y, z) coordinates.
` [ N₁,x ⋯ Nₙ,x ] [ N₁,y ⋯ Nₙ,y ] [ N₁,z ⋯ Nₙ,z ] `Obtained from
Get_dN_pg()via the inverse Jacobian.Shape:
(Ne, nPg, dim, nPe).
- Get_dN_pg(
- matrixType: MatrixType,
First derivatives of the shape functions, evaluated at Gauss points in (ξ, η, ζ).
` [ N₁,ξ ⋯ Nₙ,ξ ] [ N₁,η ⋯ Nₙ,η ] [ N₁,ζ ⋯ Nₙ,ζ ] `Shape:(nPg, dim, nPe).
- Get_ddN_e_pg(
- matrixType: MatrixType,
Second derivatives of the shape functions in physical (x, y, z) coordinates.
` [ N₁,xx ⋯ Nₙ,xx ] [ N₁,yy ⋯ Nₙ,yy ] [ N₁,zz ⋯ Nₙ,zz ] `Obtained from
Get_ddN_pg()via the (squared) inverse Jacobian.Shape:
(Ne, nPg, dim, nPe).
- Get_ddN_pg(
- matrixType: MatrixType,
Second derivatives of the shape functions, evaluated at Gauss points in (ξ, η, ζ).
` [ N₁,ξξ ⋯ Nₙ,ξξ ] [ N₁,ηη ⋯ Nₙ,ηη ] [ N₁,ζζ ⋯ Nₙ,ζζ ] `Shape:(nPg, dim, nPe).
- Get_dddN_pg(
- matrixType: MatrixType,
Third derivatives of the shape functions, evaluated at Gauss points in (ξ, η, ζ).
` [ N₁,ξξξ ⋯ Nₙ,ξξξ ] [ N₁,ηηη ⋯ Nₙ,ηηη ] [ N₁,ζζζ ⋯ Nₙ,ζζζ ] `Shape:(nPg, dim, nPe).
- Get_ddddN_pg(
- matrixType: MatrixType,
Fourth derivatives of the shape functions, evaluated at Gauss points in (ξ, η, ζ).
` [ N₁,ξξξξ ⋯ Nₙ,ξξξξ ] [ N₁,ηηηη ⋯ Nₙ,ηηηη ] [ N₁,ζζζζ ⋯ Nₙ,ζζζζ ] `Shape:(nPg, dim, nPe).
- Get_gauss(
- matrixType: MatrixType,
Returns integration points according to the matrix type.
- Get_invF_e_pg(
- matrixType: MatrixType,
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,
Returns the jacobians.
variation in size (length, area or volume) between the reference element and the actual element
- Get_leftDispPart_e_pg(
- matrixType: MatrixType,
Left half of the elastic stiffness — precomputed
wJ · Bᵀ.` leftDispPart_e_pg = wJ · Bᵀ K_e = (leftDispPart_e_pg · C · B).sum(over Gauss points) `Pulled out as a cached factor so the constitutive matrix
Cis the only operand left when assembling.Shape:
(Ne, nPg, nPe·dim, nstrain).
- Get_normals_e_pg(
- matrixType: MatrixType,
- displacementMatrix: ndarray[tuple[Any, ...], dtype[floating]] | None = None,
Returns the normals for each elements and gauss points (Ne, nPg, 3).
- Get_pointsInElem(
- coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
- elem: int,
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,
Returns the row indices used to assemble local matrices into the global matrix.
- Get_weight_pg(
- matrixType: MatrixType,
Returns integration point weights according to the matrix type.
- Get_weightedJacobian_e_pg(
- matrixType: MatrixType,
Returns the jacobian_e_pg * weight_pg.
- Integrate_e(
- func=<function _GroupElem.<lambda>>,
- matrixType=MatrixType.mass,
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,
Locates sol on elements
- Set_Tag(
- nodes: ndarray[tuple[Any, ...], dtype[int64]],
- tag: str,
Set a tag on the nodes and elements belonging to the group of elements.
- _Get_Mapping(
- coordinates_n: ndarray[tuple[Any, ...], dtype[floating]],
- elements_e: ndarray[tuple[Any, ...], dtype[int64]],
- needCoordinates=False,
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]],
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]],
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]],
Get nearby nodes.
- Parameters:
coordinates_n (_types.FloatArray) – coordinates (n, 3) array
- Returns:
nearby nodes
- Return type:
_types.IntArray
- _Get_partitioned_data() tuple[int, ndarray[tuple[Any, ...], dtype[int64]], 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,
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.
coord_e * sysCoord_e -> coordinates in element’s (i, j, k) coordinate system.
- _Init_Functions(
- order: int,
Initializes functions to be evaluated at Gauss points.
- abstractmethod _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(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]] = [],
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]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- abstractmethod _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- abstractmethod _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- abstractmethod _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- 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 center: ndarray[tuple[Any, ...], dtype[floating]][source]#
center of mass / barycenter / inertia center
- property coord: ndarray[tuple[Any, ...], dtype[floating]][source]#
this matrix contains the element group (local) coordinates (Nn, 3)
- property edges: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element edges (for FEM purposes).
- abstract property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property length_e: ndarray[tuple[Any, ...], dtype[floating]][source]#
length covered by each element
- property nodes: ndarray[tuple[Any, ...], dtype[int64]][source]#
nodes used by the element group. Node ‘n’ is on line ‘n’ in global 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.
- EasyFEA.FEM.Calc_projector( ) 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)
- EasyFEA.FEM.Det(
- mat: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
Computes det(mat)
- EasyFEA.FEM.Load_Mesh(path: str) Mesh[source]#
Loads the EasyFEA mesh from the specified pickle file.
- Parameters:
path (str) – path to the pickle file.
- Returns:
The loaded mesh.
- Return type:
- 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,
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,
np.linalg.norm() wrapper.
see https://numpy.org/doc/stable/reference/generated/numpy.linalg.norm.html
- EasyFEA.FEM.TensorProd(
- A: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- B: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- symmetric=False,
- ndim: int | None = None,
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]],
Computes trace(mat)
- EasyFEA.FEM.Transpose(
- mat: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
Computes transpose(mat)
- class EasyFEA.FEM.Elems.HEXA20(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.HEXA27(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.HEXA8(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
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 (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.POINT(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.PRISM15(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.PRISM18(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.PRISM6(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.QUAD4(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.QUAD8(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.QUAD9(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG2(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG3(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG4(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG5(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TETRA10(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.TETRA4(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- property segments: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to construct segments (for display purposes).
- class EasyFEA.FEM.Elems.TRI10(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TRI15(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TRI3(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TRI6(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElem- _N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Shape functions in (ξ, η, ζ) reference coordinates.
` [ N₁ ] [ ⋮ ] [ Nₙ ] `Shape:(nPe, 1).
- _dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
First derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξ N₁,η N₁,ζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξ Nₙ,η Nₙ,ζ ] `Shape:(nPe, dim).
- _ddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Second derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξ N₁,ηη N₁,ζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξ Nₙ,ηη Nₙ,ζζ ] `Shape:(nPe, dim).
- _dddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Third derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξ N₁,ηηη N₁,ζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξ Nₙ,ηηη Nₙ,ζζζ ] `Shape:(nPe, dim).
- _ddddN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Fourth derivatives of the shape functions in (ξ, η, ζ).
` [ N₁,ξξξξ N₁,ηηηη N₁,ζζζζ ] [ ⋮ ⋮ ⋮ ] [ Nₙ,ξξξξ Nₙ,ηηηη Nₙ,ζζζζ ] `Shape:(nPe, dim).
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Operators.Bilinear.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,
- static broadcast(
- value,
- Ne: int,
- nPg: int,
- tensor_ndim: int = 0,
Broadcast a scalar or array coefficient to a shape compatible with multiplication against an
(Ne, nPg, ...)FeArray.Returns a stride-tricked read-only view for non-scalar inputs (no data duplication). Callers must not mutate the result in place; the expected use is consumption inside expressions such as
coef * wJ_e_pg * dN_e_pg.T @ dN_e_pg, which create new arrays.tensor_ndimdeclares how many trailing axes ofvalueare tensor dims (e.g.2for a Hooke tensor(..., nstrain, nstrain)). With it set, the leading axes are checked against()/(Ne,)/(Ne, nPg)strictly — this disambiguates shapes like(Ne, n, n)from(Ne, nPg, n)whennPg == n(e.g. TRI6 in 2D, wherenPg == nstrain == 3).Accepted shapes (with default
tensor_ndim=0)#scalar (int / float / numpy scalar) → returned as
float.(Ne, nPg, ...)ndarray / FeArray → wrapped as FeArray.1-D
(Ne,)or(nPg,)→ tiled to(Ne, nPg).Any other shape broadcastable to
(Ne, nPg, ...)→ tiled with leading(Ne, nPg)dims.
- _asfearrays(
- *,
- broadcastFeArrays=False,
- dot(other, /, out=None)[source]#
Refer to
numpy.dot()for full documentation.See also
numpy.dotequivalent function
- integrate() ndarray[source]#
Integrate over the Gauss-point axis (axis 1). Returns
(Ne, ...)ndarray.
- 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
- class EasyFEA.FEM.Operators.Bilinear.MatrixType(value)[source]#
Bases:
str,EnumOrder 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
- beam_shear = 'beam_shear'#
Reduced integration order for the Timoshenko shear term (n-1 Gauss points for SEGn). Selective reduced integration is the standard cure for shear locking in Lagrange Timoshenko beams (Hughes 1987; Bathe 2006).
- mass = 'mass'#
int_Ω N • N dΩ type
- rigi = 'rigi'#
int_Ω dN • dN dΩ type
- class EasyFEA.FEM.Operators.Bilinear._EulerBernoulli(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_GroupElemEuler-Bernoulli beam element.
Kinematic assumption: plane sections remain perpendicular to the beam axis. This enforces rz = v’ (rotation = slope), so rz is NOT an independent DOF — it is derived from the transverse displacement v.
As a consequence: - shear strain γ = v’ - rz = 0 by definition (shear-rigid) - bending curvature κ = d²v/dx² (second derivative of Hermitian v) - D matrix (2D): diag([EA, EIz]) — no shear term - strain vector (2D): [du/dx, d²v/dx²] → internal forces [N, Mz]
Shape functions: - axial u: Lagrange N_i(ξ) - transverse v: cubic Hermitian Φ_i(x), Ψ_i(x) (C¹-continuous) - rotation rz: dΦ/dx, dΨ/dx (derivative of Hermitian — enslaved to v)
Shear force Ty is not a primary result; it is recovered from moment equilibrium Ty = -dMz/dx in post-processing.
Valid for slender beams (L/h ≫ 1). For stocky beams use _Timoshenko.
- Get_Hermitian_N_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Evaluates Hermitian shape functions in (x, y, z) coordinates.
[phi_i psi_i … phi_n psi_n]
(Ne, nPg, 1, nPe*2)
Phi columns are reference shape functions evaluated as-is (scalar values need no Jacobian transformation). Psi columns are scaled by L_e so that dpsi_i/dx evaluated at node i equals 1, matching the convention that the rotation DOF rz_i corresponds to v’(node i).
- Get_Hermitian_N_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#
Evaluates Hermitian shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPg, nPe*2)
- Get_Hermitian_dN_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Evaluates the first-order derivatives of Hermitian 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_Hermitian_dN_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#
Evaluates Hermitian shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPg, nPe*2)
- Get_Hermitian_ddN_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Evaluates the second-order derivatives of Hermitian 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_Hermitian_ddN_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#
Evaluates Hermitian shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξξ psi_i,ξξ … phi_n,ξξ x psi_n,ξξ]
(nPg, nPe*2)
- Get_Hermitian_dddN_e_pg() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Evaluates the third-order derivatives of Hermitian shape functions in (x, y, z) coordinates.
[phi_i,xxx psi_i,xxx … phi_n,xxx psi_n,xxx]
(Ne, nPg, 1, nPe*2)
- Get_Hermitian_dddN_pg() ndarray[tuple[Any, ...], dtype[floating]][source]#
Evaluates Hermitian shape functions third derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξξξ psi_i,ξξξ … phi_n,ξξξ x psi_n,ξξξ]
(nPg, nPe*2)
- Get_beam_B_e_pg(
- beamStructure: BeamStructure,
Euler-Bernoulli beam strain-displacement matrix B.
- Strain vector (no shear — γ = 0 by assumption):
2D: ε = [du/dx, d²v/dx²] → internal forces [N, Mz] 3D: ε = [du/dx, drx/dx, κy, κz] → internal forces [N, Mx, My, Mz]
- with κy = d²w/dx² (flex-y: Hermitian ddNvz at [w,ry] DOFs)
κz = d²v/dx² (flex-z: Hermitian ddNv at [v,rz] DOFs)
KEY difference from Timoshenko: - bending rows use ddNv (d²v/dx²), acting on BOTH displacement and rotation
DOFs simultaneously — because κ = d²v/dx² = drz/dx only under rz = v’.
there is NO shear row — B has 2 rows (2D) or 4 rows (3D).
D is 2×2 (2D) or 4×4 (3D): no kGA term.
- Get_beam_N_e_pg(
- beamStructure: BeamStructure,
Euler-Bernoulli beam shape function matrix N for the mass matrix.
Because rz = v’ (rotation enslaved to slope), the rotation row of N is filled with the Hermitian derivative dNv — NOT with independent Lagrange functions. Compare with _Timoshenko.Get_beam_N_e_pg where the rotation rows use Lagrange N.
- Get_beam_shear_B_e_pg(
- beamStructure: BeamStructure,
Shear-displacement matrix for Euler-Bernoulli shear recovery.
Ty = -dMz/dx = -EIz · d³v/dx³, Tz = -dMy/dx = -EIy · d³w/dx³. Same DOF layout and shape as Get_beam_B_e_pg, but with dddNv (third Hermitian derivative) in the bending rows. Evaluating -(D @ B_shear @ sol_e) at each Gauss point gives the exact shear for any polynomial Mz the element can represent — unlike the two-point finite difference of Mz which is only exact when Mz is linear (SEG2) or quadratic (by GP symmetry).
Returns (Ne, nPg, nStrains, dof_n * nPe) — same shape as B_e_pg.
- _Compute_P_e_pg(
- beamStructure: BeamStructure,
- abstractmethod _Hermitian_N() ndarray[tuple[Any, ...], dtype[floating]][source]#
Hermitian shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPe*2, 1)
- abstractmethod _Hermitian_dN() ndarray[tuple[Any, ...], dtype[floating]][source]#
Hermitian shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 1)
- class EasyFEA.FEM.Operators.Bilinear._Timoshenko(
- gmshId: int,
- connect: ndarray[tuple[Any, ...], dtype[int64]],
- coordinates: ndarray[tuple[Any, ...], dtype[floating]],
Bases:
_EulerBernoulliTimoshenko beam element with selective reduced integration.
Kinematic assumption: the cross-section can rotate independently of the beam slope. Rotation θ (rz in 2D, ry/rz in 3D) is an INDEPENDENT field — it is NOT constrained to equal v’. The transverse shear strain γ = v’ - θ is therefore non-zero and adds flexibility beyond bending.
Interpolation: pure Lagrange of degree (n-1) for SEGn, applied to BOTH transverse displacements (v, w) AND rotations (rx, ry, rz). Same shape functions as the axial DOF u.
Locking cure: selective reduced integration. - axial/torsion/bending terms: full Gauss (MatrixType.beam) - shear terms γ_y, γ_z : reduced Gauss (MatrixType.beam_shear, n-1 pts)
This is the standard MITC / selective-reduced-integration formulation (Hughes 1987; Bathe 2006). For tip-loaded cantilever it gives:
SEG2 (linear v): O(h²) convergence
SEG3+ (quadratic v+): machine precision at nL=1 (exact representation of the cubic bending mode)
- Strain vector:
2D: ε = [du/dx, dθz/dx, v’-θz] → [N, Mz, Ty] 3D: ε = [du/dx, drx/dx, κy, κz, γy, γz]
- with κy = -dRy/dx (sign-consistent with ry = -w’ in EB)
κz = dRz/dx γy = v’ - rz, γz = w’ + ry
D matrix (2D): diag([EA, EIz, kGA]) — 3×3, shear term added.
- Get_beam_B_e_pg(
- beamStructure: BeamStructure,
- matrixType: MatrixType = MatrixType.beam,
Timoshenko strain-displacement matrix B — pure Lagrange.
- Strain vector (shear present — γ ≠ 0):
2D: ε = [du/dx, dθz/dx, v’-θz] → [N, Mz, Ty] 3D: ε = [du/dx, drx/dx, -dRy/dx, dRz/dx, v’-rz, w’+ry]
→ [N, Mx, My, Mz, Ty, Tz]
All fields use Lagrange degree (n-1). The shear rows in particular produce equal-order locking under FULL integration, which is why the assembly (see Simulations.Beam.Construct_local_matrix_system) integrates the shear contribution at MatrixType.beam_shear (n-1 Gauss points) — selective reduced integration, the standard cure.
Passing matrixType=MatrixType.beam_shear evaluates the same B at the reduced points; passing MatrixType.beam (default) uses full points.
- Get_beam_N_e_pg(
- beamStructure: BeamStructure,
- matrixType: MatrixType = MatrixType.beam,
Timoshenko beam shape function matrix N — pure Lagrange.
Every kinematic field (u, v, w, rx, ry, rz) is interpolated with the same Lagrange polynomial of degree n-1. This is the standard Timoshenko formulation when paired with selective reduced integration of the shear term (see Get_beam_B_e_pg).
- EasyFEA.FEM.Operators.Bilinear.BeamBending(groupElem: _GroupElem, beamStructure: BeamStructure) ndarray[source]#
∫_e Bᵀ · D_bending · B dx— axial + bending (+ torsion) stiffness.Returns
(Ne, nPe·dof_n, nPe·dof_n)withdof_n = beamStructure.dof_n.Integrated at the full Gauss scheme
MatrixType.beam. For Timoshenko elements in 2D / 3D the transverse-shear rows ofDare zeroed so this term carries only the axial / bending / torsion energy; combine withBeamShear()to get the full Timoshenko stiffness, or callBeamStiffness()directly.
- EasyFEA.FEM.Operators.Bilinear.BeamMass(
- groupElem: _GroupElem,
- beamStructure: BeamStructure,
- coef: int | float | FeArray | ndarray[tuple[Any, ...], dtype[Any]] = 1.0,
∫_e coef · Nᵀ · M · N dx— beam consistent-mass matrix.Returns
(Ne, nPe·dof_n, nPe·dof_n)withdof_n = beamStructure.dof_n.Integrated at
MatrixType.beam.coefis typically the densityrhoof the simulation; may be scalar,(Ne,),(nPg,), or(Ne, nPg)— broadcast via :pymeth:`FeArray.broadcast`.
- EasyFEA.FEM.Operators.Bilinear.BeamShear(groupElem: _GroupElem, beamStructure: BeamStructure) ndarray[source]#
∫_e Bᵀ · D_shear · B dx— transverse-shear stiffness, SRI.Returns
(Ne, nPe·dof_n, nPe·dof_n)withdof_n = beamStructure.dof_n.Active only for
_TimoshenkogroupElems in 2D / 3D — returns a zero matrix otherwise. Integrated at the reduced Gauss schemeMatrixType.beam_shear(one fewer point thanMatrixType.beam) with the non-shear rows ofDzeroed. This is the selective-reduced- integration cure for Timoshenko shear locking.
- EasyFEA.FEM.Operators.Bilinear.BeamStiffness(groupElem: _GroupElem, beamStructure: BeamStructure) ndarray[source]#
Full beam stiffness
K_e = BeamBending + BeamShear.Returns
(Ne, nPe·dof_n, nPe·dof_n)withdof_n = beamStructure.dof_n.Euler-Bernoulli (or 1-D Timoshenko): reduces to
∫ Bᵀ D B dxatMatrixType.beam. Timoshenko in 2D / 3D: splits into bending (full Gauss) + shear (reduced Gauss) viaBeamBending()andBeamShear()to defeat shear locking.
- EasyFEA.FEM.Operators.Bilinear.GradUGradV(
- groupElem: _GroupElem,
- coef: int | float | FeArray | ndarray[tuple[Any, ...], dtype[Any]] = 1.0,
- matrixType: MatrixType = MatrixType.rigi,
∫_Ω coef · ∇u · ∇v dΩ— returns(Ne, nPe, nPe).coefmay be scalar,(Ne,),(nPg,), or(Ne, nPg); broadcast via :pymeth:`FeArray.broadcast` (stride view, no copy).
- EasyFEA.FEM.Operators.Bilinear.GradU_A_GradV(
- groupElem: _GroupElem,
- A: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- coef: int | float | FeArray | ndarray[tuple[Any, ...], dtype[Any]] = 1.0,
- matrixType: MatrixType = MatrixType.rigi,
∫_Ω coef · ∇u · A · ∇v dΩ— anisotropic diffusion. Returns(Ne, nPe, nPe).Ais the diffusion tensor with trailing two dims(dim, dim); the leading dims may be empty,(Ne,), or(Ne, nPg).coefis a scalar weight: scalar,(Ne,),(nPg,), or(Ne, nPg). Both broadcast via :pymeth:`FeArray.broadcast`. For the isotropic form∫ coef · ∇u · ∇v dΩ, useGradUGradV().
- EasyFEA.FEM.Operators.Bilinear.LinearizedElasticity(
- groupElem: _GroupElem,
- C: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- matrixType: MatrixType = MatrixType.rigi,
∫_Ω ε(u) : C : ε(v) dΩ— small-strain elastic stiffness.Returns
(Ne, nPe·dim, nPe·dim).Cis the Hooke tensor in Kelvin–Mandel notation with trailing two dims(nstrain, nstrain); the leading dims may be empty (homogeneous),(Ne,)(per-element), or(Ne, nPg)(per-element per-gauss-point).
- EasyFEA.FEM.Operators.Bilinear.MassAlongNormal(
- groupElem: _GroupElem,
- coef: int | float | FeArray | ndarray[tuple[Any, ...], dtype[Any]] = 1.0,
- matrixType: MatrixType = MatrixType.mass,
∫_Γ coef · (u · n̂)(v · n̂) dΓ— mass projected onto the surface normal.Returns
(Ne, nPe·3, nPe·3)in(xi, yi, zi, ..., xn, yn, zn)order.The per-Gauss-point integrand is
Nᵀ · (n̂ ⊗ n̂) · Nwheren̂is the reference-configuration unit normal supplied by_GroupElem.Get_normals_e_pg()andNis the block-diagonal shape-function matrix from_GroupElem.Get_N_pg_rep(). Typically used to penalise / enforce the normal component on a surface (Robin-stylek · u·n̂ · v·n̂boundary conditions).Restricted to a 2D surface group in a 3D mesh (
groupElem.dim == 2,groupElem.inDim == 3).coefmay be scalar,(Ne,),(nPg,), or(Ne, nPg); broadcast via :pymeth:`FeArray.broadcast` (stride view, no copy).
- EasyFEA.FEM.Operators.Bilinear.TensorProd(
- A: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- B: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- symmetric=False,
- ndim: int | None = None,
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.Operators.Bilinear.UV(
- groupElem: _GroupElem,
- coef: int | float | FeArray | ndarray[tuple[Any, ...], dtype[Any]] = 1.0,
- dof_n: int = 1,
- matrixType: MatrixType = MatrixType.mass,
∫_Ω coef · u · v dΩ— returns(Ne, nPe·dof_n, nPe·dof_n).dof_n=1for scalar fields (thermal, phase-field);dof_n=dimfor vector fields (elastic dynamics).coefmay be scalar,(Ne,),(nPg,), or(Ne, nPg); broadcast via :pymeth:`FeArray.broadcast` (stride view, no copy).
- class EasyFEA.FEM.Operators.Linear.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,
- static broadcast(
- value,
- Ne: int,
- nPg: int,
- tensor_ndim: int = 0,
Broadcast a scalar or array coefficient to a shape compatible with multiplication against an
(Ne, nPg, ...)FeArray.Returns a stride-tricked read-only view for non-scalar inputs (no data duplication). Callers must not mutate the result in place; the expected use is consumption inside expressions such as
coef * wJ_e_pg * dN_e_pg.T @ dN_e_pg, which create new arrays.tensor_ndimdeclares how many trailing axes ofvalueare tensor dims (e.g.2for a Hooke tensor(..., nstrain, nstrain)). With it set, the leading axes are checked against()/(Ne,)/(Ne, nPg)strictly — this disambiguates shapes like(Ne, n, n)from(Ne, nPg, n)whennPg == n(e.g. TRI6 in 2D, wherenPg == nstrain == 3).Accepted shapes (with default
tensor_ndim=0)#scalar (int / float / numpy scalar) → returned as
float.(Ne, nPg, ...)ndarray / FeArray → wrapped as FeArray.1-D
(Ne,)or(nPg,)→ tiled to(Ne, nPg).Any other shape broadcastable to
(Ne, nPg, ...)→ tiled with leading(Ne, nPg)dims.
- _asfearrays(
- *,
- broadcastFeArrays=False,
- dot(other, /, out=None)[source]#
Refer to
numpy.dot()for full documentation.See also
numpy.dotequivalent function
- integrate() ndarray[source]#
Integrate over the Gauss-point axis (axis 1). Returns
(Ne, ...)ndarray.
- 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
- class EasyFEA.FEM.Operators.Linear.MatrixType(value)[source]#
Bases:
str,EnumOrder 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
- beam_shear = 'beam_shear'#
Reduced integration order for the Timoshenko shear term (n-1 Gauss points for SEGn). Selective reduced integration is the standard cure for shear locking in Lagrange Timoshenko beams (Hughes 1987; Bathe 2006).
- mass = 'mass'#
int_Ω N • N dΩ type
- rigi = 'rigi'#
int_Ω dN • dN dΩ type
- EasyFEA.FEM.Operators.Linear.V(
- groupElem: _GroupElem,
- f: int | float | FeArray | ndarray[tuple[Any, ...], dtype[Any]] = 1.0,
- dof_n: int = 1,
- matrixType: MatrixType = MatrixType.mass,
∫_Ω f · v dΩ— returns(Ne, nPe·dof_n).dof_n=1for scalar fields (thermal, phase-field);dof_n=dimfor vector fields (elastic dynamics).fmay be scalar,(Ne,),(nPg,), or(Ne, nPg); broadcast via :pymeth:`FeArray.broadcast` (stride view, no copy).
- class EasyFEA.FEM.Operators.NonLinear.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,
- static broadcast(
- value,
- Ne: int,
- nPg: int,
- tensor_ndim: int = 0,
Broadcast a scalar or array coefficient to a shape compatible with multiplication against an
(Ne, nPg, ...)FeArray.Returns a stride-tricked read-only view for non-scalar inputs (no data duplication). Callers must not mutate the result in place; the expected use is consumption inside expressions such as
coef * wJ_e_pg * dN_e_pg.T @ dN_e_pg, which create new arrays.tensor_ndimdeclares how many trailing axes ofvalueare tensor dims (e.g.2for a Hooke tensor(..., nstrain, nstrain)). With it set, the leading axes are checked against()/(Ne,)/(Ne, nPg)strictly — this disambiguates shapes like(Ne, n, n)from(Ne, nPg, n)whennPg == n(e.g. TRI6 in 2D, wherenPg == nstrain == 3).Accepted shapes (with default
tensor_ndim=0)#scalar (int / float / numpy scalar) → returned as
float.(Ne, nPg, ...)ndarray / FeArray → wrapped as FeArray.1-D
(Ne,)or(nPg,)→ tiled to(Ne, nPg).Any other shape broadcastable to
(Ne, nPg, ...)→ tiled with leading(Ne, nPg)dims.
- _asfearrays(
- *,
- broadcastFeArrays=False,
- dot(other, /, out=None)[source]#
Refer to
numpy.dot()for full documentation.See also
numpy.dotequivalent function
- integrate() ndarray[source]#
Integrate over the Gauss-point axis (axis 1). Returns
(Ne, ...)ndarray.
- 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
- class EasyFEA.FEM.Operators.NonLinear.MatrixType(value)[source]#
Bases:
str,EnumOrder 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
- beam_shear = 'beam_shear'#
Reduced integration order for the Timoshenko shear term (n-1 Gauss points for SEGn). Selective reduced integration is the standard cure for shear locking in Lagrange Timoshenko beams (Hughes 1987; Bathe 2006).
- mass = 'mass'#
int_Ω N • N dΩ type
- rigi = 'rigi'#
int_Ω dN • dN dΩ type
- EasyFEA.FEM.Operators.NonLinear.FollowingPressure(
- u: ndarray,
- groupElem: _GroupElem,
- pressure: float | ndarray,
- elements: ndarray | None = None,
- matrixType: MatrixType = MatrixType.rigi,
Follower-pressure contribution on a 2D surface group in a 3D mesh.
The load tracks the deformed normal
n = ∂x/∂r × ∂x/∂swithx = X + u, so its Jacobian feeds a non-symmetric tangent.Returned
(K_e, R_e)are contributions to globalKandR(u) = R_internal − F_follower— same convention as PK2:` K_e = -∂F_follower/∂u → slot K R_e = -F_follower(u) → slot F as -R_e (= +F_follower in b) `Outside
elementsthe returned arrays are exact zero so the surface connectivity can scatter(Ne_surf, ...)uniformly.
- EasyFEA.FEM.Operators.NonLinear.KelvinVoigtDamping(
- material: _HyperElastic,
- state: HyperElasticState,
- velocity: ndarray,
Kelvin–Voigt viscous element contributions
(C_e, Kgeo_e)for the large-strain viscous forceF_visco(u) = C(u)·v, withΣ_visco = η·Ė(Green-Lagrange strain rate ofvelocity) andB = De(u)·grad.C_e = thickness · η · ∫ Bᵀ B dΩ— the damping matrix; the simulation puts it in slot 2 of(K, C, M, F)(residualb -= C @ v_t, time-scheme history, and thecoefC·Ctangent).Kgeo_e— the configuration tangent∂(C·v)/∂uat fixed velocity (geometric stiffening fromΣ_viscoplus the∂Ė/∂uterm); the simulation adds it toK_eso it ridescoefK.
- Parameters:
material – Hyperelastic constitutive law — supplies the viscosity
eta.state – Hyperelastic state — owns the mesh and the current displacement.
velocity – Velocity field (same
(xi,yi,zi,...)layout as the displacement), orNonefor a quasi-static evaluation.``(None (Returns) –
``(Ne (None)`` when material.eta == 0 or velocity is None. Both matrices are) –
nPe·dim
``(xi (nPe·dim)`` reordered to) –
yi
zi
...
xn
yn
zn)``.
- EasyFEA.FEM.Operators.NonLinear.Project_vector_to_matrix(
- vector: ndarray[tuple[Any, ...], dtype[floating]],
- coef=np.float64(1.4142135623730951),
- EasyFEA.FEM.Operators.NonLinear.SecondPiolaKirchhoffStressTensor(
- material: _HyperElastic,
- state: HyperElasticState,
Tangent and residual for a hyperelastic constitutive law.
Returns
(K_e, R_e)in(xi,yi,zi,...,xn,yn,zn).The operator pulls
De_e_pgfromstate.Compute_De()— kinematic operator,dWde_e_pgfrommaterial.Compute_dWde(state)— PK2 in Kelvin-Mandel vector form (any active-stress / viscous fold-in is the material’s responsibility),d2Wde_e_pgfrommaterial.Compute_d2Wde(state)— consistent tangent in Kelvin-Mandel matrix form,
and assembles
``` B_e_pg = De · grad (strain-displacement) Sig_e_pg = block(P(dWde_e_pg)) (geometric tangent kernel)
K_e = ∫ Bᵀ · d2Wde · B dΩ + ∫ gradᵀ · Sig · grad dΩ R_e = ∫ Bᵀ · dWde dΩ ```
where
P(·)is the Kelvin-Mandel vector → symmetric matrix projection.- Parameters:
material – Hyperelastic constitutive law — supplies
Compute_dWde(state)andCompute_d2Wde(state).state – Hyperelastic state — owns the mesh and the current displacement.
- Returns:
A_e (ndarray of shape
(Ne, nPe·dim, nPe·dim)) – Consistent tangent — sum of the linear (material) and nonlinear (geometric) pieces.r_e (ndarray of shape
(Ne, nPe·dim)) – Internal residual force.