FEM#
The EasyFEA.FEM module in EasyFEA provides essential tools for creating and managing finite element meshes, which are crucial for numerical simulations using the Finite Element Method (FEM).
What is a mesh in EasyFEA?#
A Mesh object in EasyFEA represents a collection of ElemType used to define the geometry and structure for finite element analysis. It contains multiple _GroupElem instances, which are groups of ElemType that collectively define the spatial discretization of the domain for numerical simulations.
For example, a HEXA8 mesh includes the following element types:
All implemented element types, along with their corresponding shape functions and derivatives, are defined in the EasyFEA.FEM.Elems module.
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 (e.g. .stp) part to create a mesh.Mesh_Import_mesh(): Imports an existing gmsh mesh. EasyFEA is also linked to meshio and can be used through theMedit_to_EasyFEA(): Imports medit mesh.Gmsh_to_EasyFEA(): Imports gmsh mesh.PyVista_to_EasyFEA(): Imports pyvista mesh (UnstructuredGrid or MultiBlock).Ensight_to_EasyFEA(): Imports ensight mesh.
Several examples are available in Meshes.
FEM API#
- class EasyFEA.FEM.BiLinearForm(form)[source]#
Bases:
_FormBilinear form.
- class EasyFEA.FEM.BoundaryCondition(problemType, nodes, dofs, unknowns, dofsValues, description)[source]#
Bases:
objectBoundary condition.
- static Get_dofs(problemType, list_Bc_Condition)[source]#
Returns the degrees of freedom for the problem type.
- Parameters:
problemType (str) – Problem type.
list_Bc_Condition (list[BoundaryCondition]) – List of boundary conditions.
- Returns:
degrees of freedom.
- Return type:
_types.IntArray
- static Get_dofs_nodes(availableUnknowns, nodes, unknowns)[source]#
Retrieves degrees of freedom (dofs) associated with the nodes.
- Parameters:
availableUnknowns (list[str]) – Available dofs as a list of strings. Must be a unique string list.
nodes (_types.IntArray) – Nodes for which dofs are calculated.
unknowns (list[str]) – unknowns.
- Returns:
degrees of freedom.
- Return type:
_types.IntArray
- static Get_nBc(problemType, list_Bc_Condition)[source]#
Returns the number of conditions for the problem type.
- Parameters:
problemType (str) – Problem type.
list_Bc_Condition (list[BoundaryCondition]) – List of boundary conditions.
- Returns:
Number of boundary conditions (nBc).
- Return type:
int
- static Get_values(problemType, list_Bc_Condition)[source]#
Returns the dofs values for problem type.
- Parameters:
problemType (str) – Problem type.
list_Bc_Condition (list[BoundaryCondition]) – List of boundary condition.
- Returns:
dofs values.
- Return type:
_types.FloatArray
- property dofs: ndarray[tuple[Any, ...], dtype[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 asfearray(array, broadcastFeArrays=False)[source]#
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- static ones(*shape, dtype=None)[source]#
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- static zeros(*shape, dtype=None)[source]#
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- _asfearrays(*, broadcastFeArrays=False)[source]#
- Return type:
list[Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]]
- dot(other)[source]#
Refer to
numpy.dot()for full documentation.See also
numpy.dotequivalent function
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- max(*args, **kwargs)[source]#
np.max() wrapper.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- min(*args, **kwargs)[source]#
np.min() wrapper.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- sum(*args, **kwargs)[source]#
np.sum() wrapper.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- 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, dof_n, matrixType=MatrixType.mass)[source]#
Bases:
objectField class.
- Evaluate_e(function, dofsValues, returnMeanValues=True)[source]#
Evaluates the given function for the provided field over the elements.
- Parameters:
function (Callable[[Field], FeArray]) – A function that takes a Field as input and returns a FeArray.
dofsValues (np.ndarray) – Array of shape (Nn * field.dof_n,) containing the degrees of freedom values.
returnMeanValues (bool, optional) – If True, returns the mean of the values at each element. Default is True.
- Return type:
ndarray
- Evaluate_n(function, dofsValues)[source]#
Evaluates the given function for the provided field over the nodes.
- Get_coords(concatenate=False)[source]#
Returns integration point coordinates (x,y,z) for each element.
- Interpolate(dofsValues)[source]#
Interpolates degrees of freedom values at each integration point for every element.
- Parameters:
dofsValues (np.ndarray) – Array of shape (Nn * dof_n,) containing the degrees of freedom values.
- Returns:
The (Ne, nPg, dof_n) finite element array.
- Return type:
- 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, matrixType)[source]#
Bases:
objectGauss quadrature.
- static Gauss_factory(elemType, matrixType)[source]#
Returns the integration points and weights based on the element and matrix type.
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[floating]]]
- static _Gauss_factory_nPg(elemType, nPg)[source]#
Returns the integration points and weights based on the element type and the number of Gauss points (nPg).
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[floating]]]
- static _Hexahedron(nPg)[source]#
available [8, 27]
order = [3, 5]
- Return type:
tuple[Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]]]
- static _Prism(nPg)[source]#
available [6, 8, 21]
order X = [3, 3, 5]
order Y & Z = [2, 3, 5]
- Return type:
tuple[Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]]]
- static _Quadrangle(nPg)[source]#
available [4, 9]
order = [1, 2]
- Return type:
tuple[Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]]]
- static _Tetrahedron(nPg)[source]#
available [1, 4, 5, 15]
order = [1, 2, 3, 5]
- Return type:
tuple[Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]],Collection[Union[int,float]]]
- class EasyFEA.FEM.GroupElemFactory[source]#
Bases:
object- static Create(elemType, connect, coordGlob)[source]#
Creates an element group
- Parameters:
elemType (ElemType) – element type
connect (_types.IntArray) – connection matrix storing nodes for each element (Ne, nPe)
coordGlob (_types.FloatArray) – nodes coordinates
- Returns:
the element group
- Return type:
GroupeElem
- static Get_ElemInFos(gmshId)[source]#
return elemType, nPe, dim, order, Nvertex, Nedge, Nface, Nvolume
associated with the gmsh id.
- Return type:
tuple[ElemType,int,int,int,int,int,int,int]
- static _Create(gmshId, connect, coordGlob)[source]#
Creates an element group.
- Parameters:
gmshId (int) – id gmsh
connect (_types.IntArray) – connection matrix storing nodes for each element (Ne, nPe)
coordGlob (_types.FloatArray) – nodes coordinates
- Returns:
the element group
- Return type:
GroupeElem
- static _Get_2d_element_types(elemType)[source]#
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:
gmshId
- class EasyFEA.FEM.LagrangeCondition(
- problemType,
- nodes,
- dofs,
- unknowns,
- dofsValues,
- lagrangeCoefs,
- description='',
Bases:
BoundaryConditionLagrange boundary condition
- class EasyFEA.FEM.LinearForm(form)[source]#
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()[source]#
- Return type:
list[MatrixType]
- beam = 'beam'#
int_Ω ddNv • ddNv dΩ type
- mass = 'mass'#
int_Ω N • N dΩ type
- rigi = 'rigi'#
int_Ω dN • dN dΩ type
- class EasyFEA.FEM.Mesh(dict_groupElem, verbosity=False)[source]#
Bases:
ObservableMesh class that contains several _GroupElem instances.
- Calc_regulation_projector(radius)[source]#
Returns the regulation projector matrix such that:
newU = proj * oldU
- Parameters:
radius (float) – Regularization radius for the projection operation.
- Returns:
Projection matrix of shape (Nn, Nn) that applies the regulation.
- Return type:
sp.csr_matrix
- Elements_Nodes(nodes, exclusively=True, neighborLayer=1)[source]#
Returns elements that exclusively or not use the specified nodes.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Elements_Tags(*tags)[source]#
Returns elements associated with the tag.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Evaluate_dofsValues_at_coordinates(coordinates_n, dofsValues, elements=None)[source]#
Evaluates dofsValues with shape (Nn*dof_n, ) at the specified coordinates.
- Parameters:
coordinates_n (_types.FloatArray) – coordinates that must be a (Nnodes, 3) array.
dofsValues (_types.FloatArray) – dofs values that must be a (Nn * dof_n) array.
elements (Optional[_types.IntArray], optional) – elements that may contain the specified coordinates to speed up evaluation, by default None
- Returns:
The interpolated values as a (Nnodes, dof_n) array.
- Return type:
_types.FloatArray
- Get_B_e_pg(matrixType)[source]#
Get the matrix used to calculate deformations from displacements.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
Warning
Use Kelvin Mandel Notation
[N1,x 0 … Nn,x 0
0 N1,y … 0 Nn,y
N1,y N1,x … N3,y N3,x]
(Ne, nPg, (3 or 6), nPe*dim)
- Get_DiffusePart_e_pg(matrixType)[source]#
Get the part that builds the diffusion term (scalar).
DiffusePart_e_pg = k_e_pg * jacobian_e_pg * weight_pg * dN_e_pg’ @ A @ dN_e_pg
Returns (epij) -> jacobian_e_pg * weight_pg * dN_e_pg’
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_Gradient_e_pg(u, matrixType=MatrixType.rigi)[source]#
Returns the gradient of the discretized displacement field u as a matrix
- Parameters:
u (_types.FloatArray) – discretized displacement field [ux1, uy1, uz1, …, uxN, uyN, uzN] of size Nn * dim
matrixType (MatrixType, optional) – matrix type, by default MatrixType.rigi
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]- Returns:
FeArray – grad(u) of shape (Ne, nPg, 3, 3)
dim = 1
——-
dxux 0 0
0 0 0
0 0 0
dim = 2
——-
dxux dyux 0
dxuy dyuy 0
0 0 0
dim = 3
——-
dxux dyux dzux
dxuy dyuy dzuy
dxuz dyuz dzuz
- Get_N_pg(matrixType)[source]#
Evaluates shape functions in (ξ, η, ζ) coordinates.
[N1, … , Nn]
(nPg, 1, nPe)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_N_vector_pg(matrixType)[source]#
Returns shape functions matrix in (ξ, η, ζ) coordinates
[N1 0 … Nn 0
0 N1 … 0 Nn]
(nPg, dim, npe*dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_New_meshSize_n(error_e, coef=0.5)[source]#
Returns the scalar field (at nodes) used to refine the mesh.
meshSize = (coef - 1) / error_e.max() * error_e + 1
- Parameters:
error_e (_types.FloatArray) – error evaluated on elements
coef (float, optional) – mesh size division ratio, by default 1/2
- Returns:
meshSize_n, new mesh size at nodes (Nn)
- Return type:
_types.FloatArray
- Get_Node_Values(result_e)[source]#
Get node values from element values.
The value of a node is calculated by averaging the values of the surrounding elements.
- Parameters:
mesh (Mesh) – mesh
result_e (_types.FloatArray) – element values (Ne, i)
- Returns:
nodes values (Nn, i)
- Return type:
_types.FloatArray
- Get_Paired_Nodes(corners, plot=False)[source]#
Get the paired nodes used to construct periodic boundary conditions.
- Parameters:
corners (_types.FloatArray) – Either nodes or nodes coordinates.
plot (bool, optional) – Set whether to plot the link between nodes; defaults to False.
- Returns:
Paired nodes, a 2-column matrix storing the paired nodes (n, 2).
- Return type:
_types.IntArray
- Get_Quality(criteria='aspect', nodeValues=False)[source]#
Calculates mesh quality [0, 1] (bad, good).
- Parameters:
criteria (str, optional) –
criterion used, by default ‘aspect’
”aspect”: hMin / hMax, ratio between minimum and maximum element length
”angular”: angleMin / angleMax, ratio between the minimum and maximum angle of an element
”gamma”: 2 rci/rcc, ratio between the radius of the inscribed circle and the circumscribed circle multiplied by 2. Useful for triangular elements.
”jacobian”: jMax / jMin, ratio between the maximum jacobian and the minimum jacobian. Useful for higher-order elements.
nodeValues (bool, optional) – calculates values on nodes, by default False
- Returns:
mesh quality
- Return type:
_types.FloatArray
- Get_ReactionPart_e_pg(matrixType)[source]#
Get the part that builds the reaction term (scalar).
ReactionPart_e_pg = r_e_pg * jacobian_e_pg * weight_pg * N_pg’ @ N_pg
Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’ @ N_pg
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_SourcePart_e_pg(matrixType)[source]#
Get the part that builds the source term (scalar).
SourcePart_e_pg = f_e_pg * jacobian_e_pg * weight_pg * N_pg’
Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_assembly_e(dof_n)[source]#
Returns assembly matrix for specified dof_n (Ne, nPe*dof_n)
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_columns_e(dof_n)[source]#
Returns the column indices used to assemble local matrices into the global matrix.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_connect_n_e()[source]#
Sparse matrix (Nn, Ne) of zeros and ones with ones when the node has the element such that:
values_n = connect_n_e * values_e
(Nn,1) = (Nn,Ne) * (Ne,1)
- Return type:
csr_matrix
- Get_dN_e_pg(matrixType)[source]#
Evaluates the first-order derivatives of shape functions in (x,y,z) coordinates.
[Ni,x … Nn,x
Ni,y … Nn,y]
(e, pg, dim, nPe)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_ddN_e_pg(matrixType)[source]#
Evaluates the first-order derivatives of shape functions in (x, y, z) coordinates.
[Ni,x … Nn,x
Ni,y … Nn,y
Ni,z … Nn,z]
(Ne, nPg, dim, nPe)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_jacobian_e_pg(matrixType, absoluteValues=True)[source]#
Returns the jacobians
variation in size (length, area or volume) between the reference element and the real element
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_leftDispPart(matrixType)[source]#
Get the left side of local displacement matrices.
Ku_e = jacobian_e_pg * weight_pg * B_e_pg’ @ c_e_pg @ B_e_pg
Returns (epij) -> jacobian_e_pg * weight_pg * B_e_pg’
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_list_groupElem(dim=None)[source]#
Returns the list of mesh element groups.
- Parameters:
dim (int, optional) – The dimension of elements to retrieve, by default None (uses the main mesh dimension).
- Returns:
A list of _GroupElem objects with the specified dimension.
- Return type:
list[_GroupElem]
- Get_meshSize(doMean=True)[source]#
Returns the mesh size of the mesh.
returns meshSize_e if doMean else meshSize_e_s
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_nPg(matrixType)[source]#
Returns integration points according to the matrix type.
- Return type:
int
- Get_normals(nodes=None, displacementMatrix=None)[source]#
Returns normal vectors and nodes belonging to the edge of the mesh.
returns normals, nodes.
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[int64]]]
- Get_rows_e(dof_n)[source]#
Returns the row indices used to assemble local matrices into the global matrix.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_weight_pg(matrixType)[source]#
Returns integration points according to the matrix type.
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_weightedJacobian_e_pg(matrixType)[source]#
Returns jacobian_e_pg * weight_pg.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Locates_sol_e(sol, dof_n=None, asFeArray=False)[source]#
Locates solution on elements.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Nodes_Circle(*circles, onlyOnEdge=True)[source]#
Returns the nodes in the circle(s).
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Nodes_Conditions(func)[source]#
Returns nodes that meet the specified conditions.
- Parameters:
func (function) –
Function using the x, y and z nodes coordinates and returning boolean values.
examples :
lambda x, y, z: (x < 40) & (x > 20) & (y<10)
lambda x, y, z: (x == 40) | (x == 50)
lambda x, y, z: x >= 0
- Returns:
nodes that meet the specified conditions.
- Return type:
_types.IntArray
- Nodes_Cylinder(*circles, direction=[0, 0, 1], onlyOnEdge=False)[source]#
Returns the nodes in the cylinder.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Nodes_Domain(*domains)[source]#
Returns nodes in the domain(s).
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Nodes_Line(*lines)[source]#
Returns the nodes on the line(s).
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Nodes_Point(*points)[source]#
Returns nodes on the point(s).
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Nodes_Tags(*tags)[source]#
Returns nodes associated with the tags.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Rotate(theta, center=(0, 0, 0), direction=(0, 0, 1))[source]#
Rotates the mesh coordinates around an axis.
- Parameters:
theta (float) – rotation angle [deg]
center (tuple, optional) – rotation center, by default (0,0,0)
direction (tuple, optional) – rotation direction, by default (0,0,1)
- Return type:
None
- Set_Tag(nodes, tag)[source]#
Set a tag on the nodes and elements belonging to each group of elements in the mesh.
- Symmetry(point=(0, 0, 0), n=(1, 0, 0))[source]#
Symmetrizes the mesh coordinates with respect to a specified plane.
- Parameters:
point (tuple, optional) – a point belonging to the plane, by default (0,0,0)
n (tuple, optional) – normal to the plane, by default (1,0,0)
- Return type:
None
- _Get_realistic_vector_magnitude(coef=0.1)[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
- Return type:
float
- 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 (Nn, 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=False, gmshVerbosity=False, verbosity=False)[source]#
Bases:
objectMesher class used to construct and generate the mesh via gmsh.
- static Get_Entities(points=[], lines=[], surfaces=[], volumes=[])[source]#
Get entities from from points, lines, surfaces and volumes tags
- Return type:
list[tuple[int,int]]
- static _Construct_2D_meshes(L=10, h=10, meshSize=3)[source]#
Creates 2D meshes.
- Return type:
list[Mesh]
- static _Construct_3D_meshes(L=130, h=13, b=13, meshSize=7.5, useImport3D=False)[source]#
Creates 3D meshes.
- Return type:
list[Mesh]
- Create_posFile(coord, values, folder, filename='data')[source]#
Creates of a .pos file that can be used to refine a mesh in a zone.
- Parameters:
coord (_types.FloatArray) – coordinates of values
values (_types.FloatArray) – scalar nodes values
folder (str) – save folder
filename (str, optional) – .pos file name, by default “data”.
- Returns:
Returns the path to the created .pos file
- Return type:
str
- Mesh_2D(
- contour,
- inclusions=[],
- elemType=ElemType.TRI3,
- cracks=[],
- refineGeoms=[],
- isOrganised=False,
- additionalSurfaces=[],
- additionalLines=[],
- additionalPoints=[],
- folder='',
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.
folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh
- Returns:
Created mesh
- Return type:
- Mesh_Beams(beams, elemType=ElemType.SEG2, additionalPoints=[], folder='')[source]#
Creates a beam mesh.
- Parameters:
beams (list[_Beam]) – list of Beams
elemType (ElemType, optional) – element type, by default “SEG2” [“SEG2”, “SEG3”, “SEG4”]
folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh
additionalPoints (list[Point]) – additional points that will be added to the mesh. WARNING: points must be within the domain.
- Returns:
Created mesh
- Return type:
- Mesh_Extrude(
- contour,
- inclusions=[],
- extrude=(0, 0, 1),
- layers=[],
- elemType=ElemType.TETRA4,
- cracks=[],
- refineGeoms=[],
- isOrganised=False,
- additionalSurfaces=[],
- additionalLines=[],
- additionalPoints=[],
- folder='',
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.
folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh
- Returns:
Created mesh
- Return type:
- Mesh_Import_mesh(mesh, setPhysicalGroups=False, coef=1.0)[source]#
Creates the mesh from an .msh file.
- Parameters:
mesh (str) – .msh file
setPhysicalGroups (bool, optional) – creates physical groups, by default False
coef (float, optional) – coef applied to the node coordinates, by default 1.0
- Returns:
Created mesh
- Return type:
- Mesh_Import_part(
- file,
- dim,
- meshSize=0.0,
- elemType=None,
- refineGeoms=[None],
- folder='',
Creates the mesh from .stp or .igs files.
You can only use triangles or tetrahedrons.
- Parameters:
file (str) –
.stp or .igs files.
Note that for igs files, entities cannot be recovered.
meshSize (float, optional) – mesh size, by default 0.0
elemType (ElemType, optional) – element type, by default “TRI3” or “TETRA4” depending on dim.
refineGeoms (list[Domain|Circle|str]) – geometric objects to refine the mesh
folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh
- Returns:
Created mesh
- Return type:
- Mesh_Revolve(contour, inclusions=[], axis=<EasyFEA.Geoms._line.Line object>, angle=360, layers=[30], elemType=ElemType.TETRA4, cracks=[], refineGeoms=[], isOrganised=False, additionalSurfaces=[], additionalLines=[], additionalPoints=[], folder='')[source]#
Creates a 3D mesh by rotating a surface along an axis.
- Parameters:
contour (Domain | Circle | Points | Contour) – geometry that builds the contour
inclusions (list[Domain | Circle | Points | Contour], optional) – list of hollow and filled geom objects inside the domain
axis (Line, optional) – revolution axis, by default Line(Point(), Point(0,1))
angle (float|int, optional) – revolution angle in [deg], by default 360
layers (list[int], optional) – layers in extrusion, by default [30]
elemType (ElemType, optional) – element type, by default “TETRA4” [“TETRA4”, “TETRA10”, “HEXA8”, “HEXA20”, “HEXA27”, “PRISM6”, “PRISM15”, “PRISM18”]
cracks (list[Line | Points | Contour | CircleArc]) – list of geom object used to create open or closed cracks
refineGeoms (list[Domain|Circle|str], optional) – list of geom object for mesh refinement, by default []
isOrganised (bool, optional) – mesh is organized, by default False
additionalSurfaces (list[Domain | Circle | Points | Contour]) – additional surfaces that will be added to or removed from the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). Tip: if the mesh is not well generated, you can also give the inclusions.
additionalLines (list[Union[Line,CircleArc]]) – additional lines that will be added to the surfaces created by the contour and the inclusions. (e.g Domain, Circle, Contour, Points). WARNING: lines must be within the domain.
additionalPoints (list[Point]) – additional points that will be added to the surfaces created by the contour and the inclusions. WARNING: points must be within the domain.
folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh
- Returns:
Created mesh
- Return type:
- Save_Simu(
- simu,
- results=[],
- details=False,
- edgeColor='black',
- plotMesh=True,
- showAxes=False,
- folder='',
Save the simulation in gmsh.pos format using gmsh.view
- Parameters:
simu (_Simu) – simulation
results (list[str], optional) – list of result you want to plot, by default []
details (bool, optional) – get default result values with details or not see simu.Results_nodesField_elementsField(details), by default False
edgeColor (str, optional) – color used to plot the edges, by default ‘black’
plotMesh (bool, optional) – plot the mesh, by default True
showAxes (bool, optional) – show the axes, by default False
folder (str, optional) – folder used to save .pos file, by default “”
- Return type:
None
- _Additional_Lines(dim, lines)[source]#
Adds lines to existing dim entities. WARNING: lines must be within the domain.
- _Additional_Points(dim, points, meshSize=0.0)[source]#
Adds points to existing dim entities. WARNING points must be within the domain.
- Parameters:
dim (int) – dimension (dim >= 1)
points (list[Point]) – points
meshSize (float) – meshSize
- Return type:
None
- _Additional_Surfaces(dim, surfaces, elemType, isOrganised)[source]#
Adds surfaces to existing dim entities. Tip: if the mesh is not well generated, you can also give the inclusions.
- _Cracks_SetPhysicalGroups(cracks, entities)[source]#
Creates physical groups for cracks embeded in existing gmsh entities.
returns crackLines, crackSurfaces, openPoints, openLines
- Return type:
tuple[Optional[int],Optional[int],Optional[int],Optional[int]]
- _Create_Circle(circle)[source]#
Creates a loop with a circle object.
returns loop, lines, points
- Return type:
tuple[int,list[int],list[int]]
- _Create_Contour(contour)[source]#
Creates a loop with a contour object (list of Line, CircleArc, Points).
returns loop, lines, points, openLines, openPoints
- Return type:
tuple[int,list[int],list[int],list[int],list[int]]
- _Create_Domain(domain)[source]#
Creates a loop with a domain object.
returns loop, lines, points
- Return type:
tuple[int,list[int],list[int]]
- _Create_Lines(geom, p1, p2)[source]#
Creates the lines in order to construct the contour object (Line, CircleArc, Points).
returns lines
- Return type:
list[int]
- _Extrude(surfaces, extrude=(0, 0, 1), elemType=ElemType.TETRA4, layers=[])[source]#
Extrudes gmsh surfaces and returns extruded entities.
- Parameters:
surfaces (list[int]) – gmsh surfaces
extrude (_types.Coords, optional) – extrusion vector, by default [0,0,1]
elemType (ElemType, optional) – element type used, by default “TETRA4”
layers (list[int], optional) – layers in the extrusion, by default []
- Return type:
list[tuple[int,int]]
- _Link_Contours(contour1, contour2, elemType, nLayers=0, numElems=[])[source]#
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,
- points1,
- lines2,
- points2,
- linkingLines,
- elemType,
- nLayers=0,
- listPoints=[],
Creates linking surfaces based on lines and linkingLines.
return linkingSurfaces.
- _Loop_From_Geom(geom)[source]#
Creates a loop based on the geometric object.
returns loop, lines, points
- Return type:
tuple[int,list[int],list[int]]
- _Mesh_Generate(
- dim,
- elemType,
- crackLines=None,
- crackSurfaces=None,
- openPoints=None,
- openLines=None,
- folder='',
- filename='mesh',
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
folder (str, optional) – default mesh.msh folder, by default “” does not save the mesh
filename (str, optional) – mesh saving file as filename.msh, by default mesh
- Return type:
None
- _Revolve(surfaces, axis, angle=360.0, elemType=ElemType.TETRA4, layers=[30])[source]#
Revolves gmsh surfaces and returns revolved entities.
- _Set_PhysicalGroups(
- setPoints=True,
- setLines=True,
- setSurfaces=True,
- setVolumes=True,
Creates physical groups based on created entities.
- Return type:
None
- _Set_mesh_algorithm(elemType)[source]#
Sets the mesh algorithm.
- Return type:
None
- 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)[source]#
Creates a gmsh spline from points.
returns spline, points
- Return type:
tuple[int,list[int]]
- _Surface_From_Loops(loops)[source]#
Creates a gmsh surface with a gmsh loop (must be a plane surface).
returns surface
- Return type:
int
- _Surfaces(contour, inclusions=[], elemType=ElemType.TRI3, isOrganised=False)[source]#
Creates gmsh surfaces.
They must be plane surfaces otherwise you must use ‘factory.addSurfaceFilling’ function.
returns surfaces, lines, points
- Parameters:
contour (Domain | Circle | Points | Contour) – the object that creates the surface area
inclusions (list[Domain | Circle | Points | Contour]) –
hollow or filled objects contained in the contour surface.
CAUTION: all inclusions must be contained within the contour and must not intersect.
elemType (ElemType, optional) – element type used, by default TRI3
isOrganised (bool, optional) – mesh is organized, by default False
- Return type:
tuple[list[int],list[int],list[int]]
- _Surfaces_Organize(surfaces, elemType, isOrganised=False, numElems=[])[source]#
Organizes surfaces.
- Parameters:
surfaces (list[int]) – list of gmsh surfaces
elemType (ElemType) – element type
isOrganised (bool, optional) – the mesh is organized, by default False
numElems (_types.Numbers, optional) – number of elements per line, by default []
- Return type:
None
- class EasyFEA.FEM._Form(form)[source]#
Bases:
ABCForm class from which BilinearForm and LinearForm are derived.
- class EasyFEA.FEM._GroupElem(gmshId, connect, coordGlob)[source]#
Bases:
ABCThe _GroupElem base class, from which all element types inherit.
- static _Eval_Functions(functions, gaussPoints)[source]#
Evaluates functions at coordinates.
Use this function to evaluate shape functions.
- Parameters:
functions (_types.FloatArray) – functions to evaluate, (nPe, dim)
gaussPoints (_types.FloatArray) –
gauss coordinates where functions will be evaluated (nPg, dim).
Be careful dim must be consistent with function arguments
- Returns:
Evaluated functions (nPg, dim, nPe)
- Return type:
_types.FloatArray
- static _Get_assembly_e(connect, dof_n)[source]#
Get the assembly matrix for the specified dof_n (Ne, nPe*dof_n)
- Parameters:
connect (_types.IntArray) – connectivity matrix (Ne, nPe)
dof_n (int) – degree of freedom per node
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_B_e_pg(matrixType)[source]#
Get the matrix used to calculate deformations from displacements.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
Warning
Use Kelvin Mandel Notation
[N1,x 0 … Nn,x 0
0 N1,y … 0 Nn,y
N1,y N1,x … N3,y N3,x]
(Ne, nPg, (3 or 6), nPe*dim)
- Get_DiffusePart_e_pg(matrixType)[source]#
Get the part that builds the diffusion term (scalar).
DiffusePart_e_pg = k_e_pg * jacobian_e_pg * weight_pg * dN_e_pg’ @ A @ dN_e_pg
Returns (epij) -> jacobian_e_pg * weight_pg * dN_e_pg’
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_Elements_Nodes(nodes, exclusively=True)[source]#
Returns elements that exclusively or not use the specified nodes.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_Elements_Tag(tag)[source]#
Returns elements associated with the tag.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_EulerBernoulli_B_e_pg(beamStructure)[source]#
Get Euler-Bernoulli beam shape functions derivatives
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_EulerBernoulli_N_e_pg()[source]#
Evaluates Euler-Bernoulli beam shape functions in (x, y, z) coordinates.
[phi_i psi_i … phi_n psi_n]
(Ne, nPg, 1, nPe*2)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_EulerBernoulli_N_e_pg_for_beam(beamStructure)[source]#
Euler-Bernoulli beam shape functions.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_EulerBernoulli_N_pg()[source]#
Evaluates Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPg, nPe*2)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_EulerBernoulli_dN_e_pg()[source]#
Evaluates the first-order derivatives of Euler-Bernoulli beam shape functions in (x, y, z) coordinates.
[phi_i,x psi_i,x … phi_n,x psi_n,x]
(Ne, nPg, 1, nPe*2)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_EulerBernoulli_dN_pg()[source]#
Evaluates Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPg, nPe*2)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_EulerBernoulli_ddN_e_pg()[source]#
Evaluates the second-order derivatives of Euler-Bernoulli beam shape functions in (x, y, z) coordinates.
[phi_i,xx psi_i,xx … phi_n,xx psi_n,xx]
(Ne, nPg, 1, nPe*2)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_EulerBernoulli_ddN_pg()[source]#
Evaluates Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ x psi_n,ξ]
(nPg, nPe*2)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_F_e_pg(matrixType)[source]#
Returns the transposed Jacobian matrix.
This matrix describes the transformation of the (ξ, η, ζ) axes from the reference element to the (x, y, z) coordinate system of the actual element.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_GaussCoordinates_e_pg(matrixType, elements=array([], dtype=float64))[source]#
Returns integration point coordinates for each element (Ne, nPg, 3) in the (x, y, z) coordinates.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_Gradient_e_pg(u, matrixType=MatrixType.rigi)[source]#
Returns the gradient of the discretized displacement field u as a matrix
- Parameters:
u (_types.FloatArray) – discretized displacement field [ux1, uy1, uz1, …, uxN, uyN, uzN] of size Nn * dim
matrixType (MatrixType, optional) – matrix type, by default MatrixType.rigi
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]- Returns:
FeArray – grad(u) of shape (Ne, nPg, 3, 3)
dim = 1
——-
dxux 0 0
0 0 0
0 0 0
dim = 2
——-
dxux dyux 0
dxuy dyuy 0
0 0 0
dim = 3
——-
dxux dyux dzux
dxuy dyuy dzuy
dxuz dyuz dzuz
- abstractmethod Get_Local_Coords()[source]#
Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_Mapping(coordinates_n, elements_e=None, needCoordinates=False)[source]#
Locates coordinates within elements.
- Return type:
tuple[ndarray[tuple[Any,...],dtype[int64]],ndarray[tuple[Any,...],dtype[int64]],ndarray[tuple[Any,...],dtype[int64]],Optional[ndarray[tuple[Any,...],dtype[floating]]]]- Returns:
detectedNodes : The nodes that have been identified within the detected elements with shape=(Nn,).
detectedElements_e : The elements in which the nodes have been detected with shape=(Ne,).
- connect_e_nThe connectivity matrix that includes the nodes identified in each element with shape=(Ne, ?).
The “?” indicates that the array may have varying dimensions along axis=1.
- coordInElem_nThe coordinates of the identified nodes, expressed in the reference element’s (ξ, η, ζ) coordinate system.
This is applicable only if needCoordinates is set to True.
- Get_N_pg(matrixType)[source]#
Evaluates shape functions in (ξ, η, ζ) coordinates.
[N1, … , Nn]
(nPg, 1, nPe)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_N_pg_rep(matrixType, repeat=1)[source]#
Repeats shape functions in the (ξ, η, ζ) coordinates.
- Parameters:
matrixType (MatrixType) – matrix type
repeat (int, optional) – number of repetitions, by default 1
Returns
-------
(nPg (• Scalar shape functions) –
[Ni 0 … Nn 0
0 Ni … 0 Nn]
rep=2 –
[Ni 0 … Nn 0
0 Ni … 0 Nn]
rep=2*dim) –
[Ni 0 … Nn 0
0 Ni … 0 Nn]
(nPg – [Ni … Nn]
rep=1 – [Ni … Nn]
nPe) – [Ni … Nn]
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_Nodes_Circle(circle, onlyOnEdge=False)[source]#
Returns nodes in the circle.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_Nodes_Conditions(func)[source]#
Returns nodes that meet the specified conditions.
- Parameters:
func (
Callable) –Function using x, y and z nodes coordinates and returning boolean values.
examples :
lambda x, y, z: (x < 40) & (x > 20) & (y<10)
lambda x, y, z: (x == 40) | (x == 50)
lambda x, y, z: x >= 0
- Returns:
nodes that meet conditions
- Return type:
_types.IntArray
- Get_Nodes_Cylinder(circle, direction=[0, 0, 1], onlyOnEdge=False)[source]#
Returns nodes in the cylinder.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_Nodes_Domain(domain)[source]#
Returns nodes in the domain.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_Nodes_Line(line)[source]#
Returns nodes on the line.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_Nodes_Point(point)[source]#
Returns nodes on the point.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_Nodes_Tag(tag)[source]#
Returns node associated with the tag.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_ReactionPart_e_pg(matrixType)[source]#
Get the part that builds the reaction term (scalar).
ReactionPart_e_pg = r_e_pg * jacobian_e_pg * weight_pg * N_pg’ @ N_pg
Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’ @ N_pg
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_SourcePart_e_pg(matrixType)[source]#
Get the part that builds the source term (scalar).
SourcePart_e_pg = f_e_pg * jacobian_e_pg * weight_pg * N_pg’
Returns (epij) -> jacobian_e_pg * weight_pg * N_pg’
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_assembly_e(dof_n)[source]#
Get the assembly matrix for the specified dof_n (Ne, nPe*dof_n)
- Parameters:
dof_n (int) – degree of freedom per node
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_columns_e(dof_n)[source]#
Returns the column indices used to assemble local matrices into the global matrix.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_connect_n_e()[source]#
Sparse matrix (Nn, Ne) of zeros and ones with ones when the node has the element such that: values_n = connect_n_e * values_e
(Nn,1) = (Nn,Ne) * (Ne,1)
- Return type:
csr_matrix
- Get_dN_e_pg(matrixType)[source]#
Evaluates the first-order derivatives of shape functions in (x, y, z) coordinates.
[Ni,x … Nn,x
Ni,y … Nn,y
Ni,z … Nn,z]
(Ne, nPg, dim, nPe)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_dN_pg(matrixType)[source]#
Evaluates shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ … Nn,ξ
Ni,η … Nn,η
Ni,ζ … Nn,ζ
(nPg, dim, nPe)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_ddN_e_pg(matrixType)[source]#
Evaluates the second-order derivatives of shape functions in (x, y, z) coordinates.
[Ni,x2 … Nn,x2
Ni,y2 … Nn,y2
Ni,z2 … Nn,z2]
(Ne, nPg, dim, nPe)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_ddN_pg(matrixType)[source]#
Evaluates shape functions second derivatives in the (ξ, η, ζ) coordinates.
[Ni,ξ2 … Nn,ξ2
Ni,η2 … Nn,η2
Ni,ζ2 … Nn,ζ2]
(nPg, dim, nPe)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_dddN_pg(matrixType)[source]#
Evaluates shape functions third derivatives in the (ξ, η, ζ) coordinates.
[Ni,ξ3 … Nn,ξ3
Ni,η3 … Nn,η3
Ni,ζ3 … Nn,ζ3]
(nPg, dim, nPe)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_ddddN_pg(matrixType)[source]#
Evaluates shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
[Ni,ξ4 … Nn,ξ4
Ni,η4 … Nn,η4
Ni,ζ4 … Nn,ζ4]
(pg, dim, nPe)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_gauss(matrixType)[source]#
Returns integration points according to the matrix type.
- Return type:
- Get_invF_e_pg(matrixType)[source]#
Returns the inverse of the transposed Jacobian matrix.
Used to obtain the derivative of the dN_e_pg shape functions in the actual element dN_e_pg = invF_e_pg • dN_pg
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_jacobian_e_pg(matrixType, absoluteValues=True)[source]#
Returns the jacobians.
variation in size (length, area or volume) between the reference element and the actual element
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_leftDispPart(matrixType)[source]#
Get the left side of local displacement matrices.
Ku_e = jacobian_e_pg * weight_pg * B_e_pg’ @ c_e_pg @ B_e_pg
Returns (epij) -> jacobian_e_pg * weight_pg * B_e_pg’
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_normals_e_pg(matrixType, displacementMatrix=None)[source]#
Returns the normals for each elements and gauss points (Ne, nPg, 3).
- Return type:
- Get_pointsInElem(coordinates_n, elem)[source]#
Returns the indexes of the coordinates contained in the element.
- Parameters:
coordinates_n (_types.FloatArray) – coordinates (n, 3)
elem (int) – element
- Returns:
indexes of coordinates contained in element
- Return type:
_types.IntArray
- Get_rows_e(dof_n)[source]#
Returns the row indices used to assemble local matrices into the global matrix.
- Return type:
ndarray[tuple[Any,...],dtype[int64]]
- Get_weight_pg(matrixType)[source]#
Returns integration point weights according to the matrix type.
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- Get_weightedJacobian_e_pg(matrixType)[source]#
Returns the jacobian_e_pg * weight_pg.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Integrate_e(func=<function _GroupElem.<lambda>>, matrixType=MatrixType.mass)[source]#
Integrates the function over elements.
- Parameters:
func (lambda) –
function that uses the x,y,z coordinates of the element’s integration points
Examples:
lambda x,y,z: 1 -> that will just integrate the element lambda x,y,z: x lambda x,y,z: x + y
lambda x,y,z: z**2
matrixType (MatrixType, optional) – matrix type, by default MatrixType.mass
- Returns:
integrated values on elements
- Return type:
_types.FloatArray
- Locates_sol_e(sol, dof_n=None, asFeArray=False)[source]#
Locates sol on elements
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Set_Tag(nodes, tag)[source]#
Set a tag on the nodes and elements belonging to the group of elements.
- _EulerBernoulli_N()[source]#
Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _EulerBernoulli_dN()[source]#
Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _EulerBernoulli_ddN()[source]#
Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 2)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _Get_Mapping(coordinates_n, elements_e, needCoordinates=False)[source]#
Locates coordinates within elements.
- Return type:
tuple[ndarray[tuple[Any,...],dtype[int64]],ndarray[tuple[Any,...],dtype[int64]],ndarray[tuple[Any,...],dtype[int64]],Optional[ndarray[tuple[Any,...],dtype[floating]]]]- Returns:
detectedNodes : The nodes that have been identified within the detected elements with shape=(Nn,).
detectedElements_e : The elements in which the nodes have been detected with shape=(Ne,).
- connect_e_nThe connectivity matrix that includes the nodes identified in each element with shape=(Ne, ?).
The “?” indicates that the array may have varying dimensions along axis=1.
- coordInElem_nThe coordinates of the identified nodes, expressed in the reference element’s (ξ, η, ζ) coordinate system.
This is applicable only if needCoordinates is set to True.
- _Get_coord_Near(coordinates_n, coordElem, dims)[source]#
Get indexes in coordinates_n that are within the coordElem’s bounds.
- Parameters:
coordinates_n (_types.FloatArray) – coordinates to check
coordElem (_types.FloatArray) – element’s bounds
dims (_types.FloatArray) – (nX, nY, nZ) = np.max(coordinates_n, 0) - np.min(coordinates_n, 0) + 1
- Returns:
indexes in element’s bounds.
- Return type:
_types.IntArray
- _Get_nearby_elements(coordinates_n)[source]#
Get nearby elements.
- Parameters:
coordinates_n (_types.FloatArray) – coordinates (n, 3) array
- Returns:
nearby elements
- Return type:
_types.IntArray
- _Get_nearby_nodes(coordinates_n)[source]#
Get nearby nodes.
- Parameters:
coordinates_n (_types.FloatArray) – coordinates (n, 3) array
- Returns:
nearby nodes
- Return type:
_types.IntArray
- _Get_partitionned_data()[source]#
Returns the paritionned data used in mpi.
(rank, elementsGlob, nodes, ghostNodes)
- Return type:
tuple[ndarray[tuple[Any,...],dtype[int64]],ndarray[tuple[Any,...],dtype[int64]],ndarray[tuple[Any,...],dtype[int64]]]
- _Get_sysCoord_e(displacementMatrix=None)[source]#
Get the basis transformation matrix (Ne, 3, 3).
[ix, jx, kx
iy, jy, ky
iz, jz, kz]
This matrix can be used to project points with (x, y, z) coordinates into the element’s (i, j, k) coordinate system.
coordo_e * sysCoord_e -> coordinates in element’s (i, j, k) coordinate system.
- _InitMatrix()[source]#
Initializes matrix dictionaries for finite element construction
- Return type:
None
- _Init_Functions(order)[source]#
Initializes functions to be evaluated at Gauss points.
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- abstractmethod _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _Set_partitionned_data(elementsGlob, nodes, rank=0)[source]#
Sets the paritionned data used in mpi.
- Parameters:
elementsGlob (_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
Remark
------
array. (Ghost nodes will be computed using the given (non-ghost) nodes)
- Return type:
None
- abstractmethod _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- abstractmethod _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- abstractmethod _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- abstractmethod _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property _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 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 coordinates (Nn, 3)
- property coordGlob: ndarray[tuple[Any, ...], dtype[floating]][source]#
this matrix contains all the mesh coordinates (mesh.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 coordGlob
- 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(oldMesh, newMesh)[source]#
Get the matrix used to project the solution from the old mesh to the new mesh such that:
newU = proj * oldU
(newNn) = (newNn x oldNn) (oldNn)
- EasyFEA.FEM.Det(mat)[source]#
Computes det(mat)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- EasyFEA.FEM.Mesh_Optim(
- DoMesh,
- folder,
- criteria='aspect',
- quality=0.8,
- ratio=0.7,
- iterMax=20,
- coef=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, **kwargs)[source]#
np.linalg.norm() wrapper.
see https://numpy.org/doc/stable/reference/generated/numpy.linalg.norm.html
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- EasyFEA.FEM.TensorProd(A, B, symmetric=False, ndim=None)[source]#
Computes tensor product.
- Parameters:
A (FeArray.FeArrayALike) – array A
B (FeArray.FeArrayALike) – array B
symmetric (bool, optional) – do symmetric product, by default False
ndim (int, optional) – ndim=1 -> vect or ndim=2 -> matrix, by default None
- Returns:
the calculated tensor product
- Return type:
FeArray.FeArrayALike
- EasyFEA.FEM.Trace(mat)[source]#
Computes trace(mat)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- EasyFEA.FEM.Transpose(mat)[source]#
Computes transpose(mat)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- class EasyFEA.FEM.Elems.HEXA20(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- Get_Local_Coords()[source]#
Get local ξ, η, ζ coordinates as a (nPe, dim) numpy array
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.PRISM15(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.QUAD8(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.QUAD9(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG2(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _EulerBernoulli_N()[source]#
Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_dN()[source]#
Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_ddN()[source]#
Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 2)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG3(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _EulerBernoulli_N()[source]#
Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_dN()[source]#
Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_ddN()[source]#
Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 2)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG4(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _EulerBernoulli_N()[source]#
Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_dN()[source]#
Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_ddN()[source]#
Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 2)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.SEG5(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _EulerBernoulli_N()[source]#
Euler-Bernoulli beam shape functions in the (ξ, η, ζ) coordinates.
[phi_i psi_i … phi_n psi_n]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_dN()[source]#
Euler-Bernoulli beam shape functions first derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 1)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _EulerBernoulli_ddN()[source]#
Euler-Bernoulli beam shape functions second derivatives in the (ξ, η, ζ) coordinates.
[phi_i,ξ psi_i,ξ … phi_n,ξ psi_n,ξ]
(nPe*2, 2)
- Return type:
ndarray[tuple[Any,...],dtype[Any]]
- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TETRA10(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[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, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TRI15(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TRI3(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).
- class EasyFEA.FEM.Elems.TRI6(gmshId, connect, coordGlob)[source]#
Bases:
_GroupElem- _N()[source]#
Shape functions in (ξ, η, ζ) coordinates.
N1
⋮
Nn
(nPe, 1)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dN()[source]#
Shape functions first derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ Ni,η Ni,ζ
⋮
Nn,ξ Nn,η Nn,ζ
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddN()[source]#
Shape functions second derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ2 Ni,η2 Ni,ζ2
⋮
Nn,ξ2 Nn,η2 Nn,ζ2
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _dddN()[source]#
Shape functions third derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ3 Ni,η3 Ni,ζ3
⋮
Nn,ξ3 Nn,η3 Nn,ζ3
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- _ddddN()[source]#
Shape functions fourth derivatives in the (ξ, η, ζ) coordinates.
Ni,ξ4 Ni,η4 Ni,ζ4
⋮
Nn,ξ4 Nn,η4 Nn,ζ4
(nPe, dim)
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property faces: ndarray[tuple[Any, ...], dtype[int64]][source]#
array of indices used to form the element faces (for FEM purposes).