Models#
The EasyFEA.Models module provides essential tools for creating and managing models.
These models are used to build _Simu instances and mainly contain material parameters.
In the simulation workflow, Models is the third step: once the mesh exists, a model encapsulates the physics and material constants (Young’s modulus, thermal conductivity, fracture toughness, …) before being passed to the simulation constructor.
With this module, you can construct:
Linear elastic materials, such as
Isotropic,TransverselyIsotropic,Orthotropic, andAnisotropic, inModels.Elastic.
Nonlinear hyperelastic materials, such as
NeoHookean,CiarletGeymonat,MooneyRivlin,SaintVenantKirchhoff, andHolzapfelOgden, inModels.HyperElastic.
Elastic beams with
Isotropic,BeamStructure, inModels.Beam.Phase-field materials with
PhaseField.Thermal materials with
Thermal.Weak forms with
WeakForms.
See also
Models API#
Module implementing constitutive laws used in simulations.
- class EasyFEA.Models.ModelType(value)[source]#
Bases:
str,EnumModel types.
- beam = 'beam'#
- damage = 'damage'#
- elastic = 'elastic'#
- hyperelastic = 'hyperelastic'#
- thermal = 'thermal'#
- weakForm = 'weakForm'#
- class EasyFEA.Models.PhaseField(
- material: _Elastic,
- split: SplitType,
- regularization: ReguType,
- Gc: float | ndarray[tuple[Any, ...], dtype[floating]],
- l0: float,
- solver=SolverType.History,
- A: ndarray[tuple[Any, ...], dtype[floating]] | None = None,
Bases:
_IModelPhaseField class.
- class SolverType(value)[source]#
Bases:
str,EnumSolver used to manage crack irreversibility.
- BoundConstrain = 'BoundConstrain'#
- History = 'History'#
- HistoryDamage = 'HistoryDamage'#
- class SplitType(value)[source]#
Bases:
str,EnumSplit Models.
- Amor = 'Amor'#
- AnisotStrain = 'AnisotStrain'#
- AnisotStrain_MP = 'AnisotStrain_MP'#
- AnisotStrain_NoCross = 'AnisotStrain_NoCross'#
- AnisotStrain_PM = 'AnisotStrain_PM'#
- AnisotStress = 'AnisotStress'#
- AnisotStress_MP = 'AnisotStress_MP'#
- AnisotStress_NoCross = 'AnisotStress_NoCross'#
- AnisotStress_PM = 'AnisotStress_PM'#
- Bourdin = 'Bourdin'#
- He = 'He'#
- Miehe = 'Miehe'#
- Stress = 'Stress'#
- Zhang = 'Zhang'#
- static Get_solvers() list[SolverType][source]#
Returns available solvers used to manage crack irreversibility
- Calc_C(
- Epsilon_e_pg: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- verif=False,
Computes the splited stifness matrices for the given strain field.
- Parameters:
Epsilon_e_pg (FeArray.FeArrayALike) – strains field (e, p, D)
- Returns:
cP_e_pg, cM_e_pg: positive and negative stifness matrices (e, p, D, D)
- Return type:
- Calc_Sigma_e_pg(
- Epsilon_e_pg: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
Computes the Stress field using the strains and the split such that:
SigmaP_e_pg = cP_e_pg * Epsilon_e_pg
SigmaM_e_pg = cM_e_pg * Epsilon_e_pg
- Parameters:
Epsilon_e_pg (FeArray.FeArrayALike) – strains field (e, p, D)
- Returns:
SigmaP_e_pg, SigmaM_e_pg: positive and negative stress fields (e, p, D)
- Return type:
- Calc_psi_e_pg(
- Epsilon_e_pg: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
Computes the elastic energy densities.
psiP_e_pg = 1/2 SigmaP_e_pg * Epsilon_e_pg
psiM_e_pg = 1/2 SigmaM_e_pg * Epsilon_e_pg
Such as :
SigmaP_e_pg = cP_e_pg * Epsilon_e_pg
SigmaM_e_pg = cM_e_pg * Epsilon_e_pg
- Get_f_e_pg(
- PsiP_e_pg: ndarray[tuple[Any, ...], dtype[floating]],
Returns source therm
- Get_g_e_pg(
- d_n: ndarray[tuple[Any, ...], dtype[floating]],
- mesh: Mesh,
- matrixType: MatrixType,
- k_res=1e-12,
Returns degradation function
- Get_r_e_pg(
- PsiP_e_pg: ndarray[tuple[Any, ...], dtype[floating]],
Returns reaction therm
- _Eigen_values_vectors_projectors(
- vector_e_pg: FeArray | ndarray[tuple[Any, ...], dtype[Any]],
- verif=False,
Computes the eigen values and eigen projectors of a second-order tensor (as a vector).
- A#
matrix characterizing the weak anisotropy in the crack surface density function
- Gc: float#
critical energy release rate
- l0: float#
half crack width
- solver: SolverType#
solver used to manage crack irreversibility
- class EasyFEA.Models.Thermal(k: float, c=0.0, thickness: float = 1.0)[source]#
Bases:
_IModelThermal class.
- c#
mass heat capacity [J K^-1 kg^-1]
- dim: int#
- k#
thermal conductivity [W m^-1]
- thickness: float#
- class EasyFEA.Models.WeakForms(
- field: Field,
- computeK: BiLinearForm,
- computeC: BiLinearForm | None = None,
- computeM: BiLinearForm | None = None,
- computeF: LinearForm | None = None,
- thickness: float = 1.0,
Bases:
_IModelClass responsible for computing the finite element matrices used in the system \(\Krm \, \mathrm{u} + \Crm \, \vrm + \Mrm \, \arm = \Frm\).
- property computeC: BiLinearForm | None[source]#
Function used to build stiffness matrix \(\Crm\).
- property computeF: LinearForm | None[source]#
Function used to build force vector \(\Frm\).
- property computeK: BiLinearForm | None[source]#
Function used to build stiffness matrix \(\Krm\).
- property computeM: BiLinearForm | None[source]#
Function used to build stiffness matrix \(\Mrm\).
- thickness: float#
- class EasyFEA.Models._IModel[source]#
Bases:
Observable,Updatable,ABCModel interface.
- EasyFEA.Models.Apply_Pmat(
- P: ndarray[tuple[Any, ...], dtype[floating]],
- M: ndarray[tuple[Any, ...], dtype[floating]],
- toGlobal=True,
Performs a basis transformation from the material’s coordinate system to the (x,y,z) coordinate system to orient the material in space.
Caution: P must be in Kelvin mandel notation
- Parameters:
P (_types.FloatArray) – P in mandel notation obtained with Get_Pmat
M (_types.FloatArray) – 3x3 or 6x6 matrix
toGlobal (bool, optional) –
sets wheter you want to get matrix in global or material coordinates, by default True
if toGlobal:
Matrix_global = P * C_material * P’
else:
Matrix_material = P’ * Matrix_global * P
- Returns:
new matrix
- Return type:
_types.FloatArray
- EasyFEA.Models.Get_Pmat(
- axis_1: ndarray[tuple[Any, ...], dtype[floating]],
- axis_2: ndarray[tuple[Any, ...], dtype[floating]],
- useMandel=True,
Constructs Pmat to pass from the material coordinates (x,y,z) to the global coordinate (X,Y,Z) such that:
if useMandel:
return [Pm]
else:
return [Ps], [Pe]
In Kelvin Mandel notation:#
Sig & Eps en [11, 22, 33, sqrt(2)*23, sqrt(2)*13, sqrt(2)*12]
[C_global] = [Pm] * [C_material] * [Pm]^T & [C_material] = [Pm]^T * [C_global] * [Pm]
[S_global] = [Pm] * [S_material] * [Pm]^T & [S_material] = [Pm]^T * [S_global] * [Pm]
Sig_global = [Pm] * Sig_material & Sig_material = [Pm]^T * Sig_global
Eps_global = [Pm] * Eps_material & Eps_material = [Pm]^T * Eps_global
In Voigt’s notation:#
Sig [S11, S22, S33, S23, S13, S12]
Eps [E11, E22, E33, 2*E23, 2*E13, 2*E12]
[C_global] = [Ps] * [C_material] * [Ps]^T & [C_material] = [Pe]^T * [C_global] * [Pe]
S_global = [Pe] * [S_material] * [Pe]^T & [S_material] = [Ps]^T * S_global * [Ps]
Sig_global = [Ps] * Sig_material & Sig_material = [Pe]^T * Sig_global
Eps_global = [Pe] * Eps_material & Eps_material = [Ps]^T * Eps_global
P matrices are orhogonal such that: inv([P]) = [P]^T
Here we use “Chevalier 1988 : Comportements élastique et viscoélastique des composites”
- EasyFEA.Models.Heterogeneous_Array(
- array: ndarray[tuple[Any, ...], dtype[floating]],
Builds a heterogeneous array.
- EasyFEA.Models.KelvinMandel_Matrix(
- dim: int,
- M: ndarray[tuple[Any, ...], dtype[floating]],
Apply Kelvin Mandel coefficient to constitutive laws.
In 2D:#
[1, 1, r2]
[1, 1, r2]
[r2, r2, 2]]
In 3D:#
[1,1,1,r2,r2,r2]
[1,1,1,r2,r2,r2]
[1,1,1,r2,r2,r2]
[r2,r2,r2,2,2,2]
[r2,r2,r2,2,2,2]
[r2,r2,r2,2,2,2]]
- EasyFEA.Models.Project_Kelvin(
- A: ndarray[tuple[Any, ...], dtype[floating]],
- orderA: int | None = None,
Projects the tensor A in Kelvin Mandel notation.
- Parameters:
A (_types.FloatArray) – tensor A (2 or 4 order tensor)
orderA (int, optional) – tensor order, by default None
- Returns:
Projected tensor
- Return type:
_types.FloatArray
- EasyFEA.Models.Project_matrix_to_vector(
- matrix: ndarray[tuple[Any, ...], dtype[floating]],
- coef=np.float64(1.4142135623730951),
- EasyFEA.Models.Project_vector_to_matrix(
- vector: ndarray[tuple[Any, ...], dtype[floating]],
- coef=np.float64(1.4142135623730951),
- EasyFEA.Models.Reshape_variable(
- variable: int | float | ndarray[tuple[Any, ...], dtype[Any]],
- Ne: int,
- nPg: int,
Resizes variable to (Ne, nPg, …) shape.
- EasyFEA.Models.Result_in_Strain_or_Stress_field(
- field_e_pg: FeArray,
- result: str,
- coef=np.float64(1.4142135623730951),
Extracts a specific result from a 2D or 3D strain or stress field.
- Parameters:
field_e_pg (_types.FloatArray) – Strain or stress field in each element and gauss points.
result (str) –
Desired result/value to extract:
2D: [xx, yy, xy, vm, Strain, Stress]
3D: [xx, yy, zz, yz, xz, xy, vm, Strain, Stress]
coef (float, optional) – Coefficient used to scale cross components in the field (e.g., xy/coef in dim=2, or yz/coef, xz/coef, xy/coef if dim=3).
- Returns:
The extracted field corresponding to the specified result.
- Return type:
_types.FloatArray
Elastic module.
- class EasyFEA.Models.Elastic.Anisotropic(
- dim: int,
- C: ndarray[tuple[Any, ...], dtype[floating]],
- useVoigtNotation: bool,
- axis1: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (1, 0, 0),
- axis2: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 1, 0),
- thickness=1.0,
Bases:
_ElasticAnisotropic Linearized Elastic material.
- Set_C(
- C: ndarray[tuple[Any, ...], dtype[floating]],
- useVoigtNotation=True,
- update_S=True,
Updates the constitutives laws by updating the C stiffness and S compliance matrices in Kelvin Mandel notation.
- Parameters:
C (_types.FloatArray) – Stifness matrix (Lamé’s law)
useVoigtNotation (bool, optional) – uses Kevin Mandel’s notation, by default True
update_S (bool, optional) – updates the compliance matrix (Hooke’s law), by default True
- Walpole_Decomposition() tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- _Behavior(
- C: ndarray[tuple[Any, ...], dtype[floating]],
- useVoigtNotation: bool,
- class EasyFEA.Models.Elastic.Isotropic(dim: int, E=210000.0, v=0.3, planeStress=True, thickness=1.0)[source]#
Bases:
_ElasticIsotropic Linearized Elastic material.
- Walpole_Decomposition() tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- _Behavior(dim: int | None = None)[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices in Kelvin Mandel notation.
In 2D:#
C -> C : Epsilon = Sigma [Sxx Syy sqrt(2)*Sxy]
S -> S : Sigma = Epsilon [Exx Eyy sqrt(2)*Exy]
In 3D:#
C -> C : Epsilon = Sigma [Sxx Syy Szz sqrt(2)*Syz sqrt(2)*Sxz sqrt(2)*Sxy]
S -> S : Sigma = Epsilon [Exx Eyy Ezz sqrt(2)*Eyz sqrt(2)*Exz sqrt(2)*Exy]
- _Update() None[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices. in Kelvin Mandel notation
- E: float#
Young’s modulus
- v: float#
Poisson’s ratio (-1<v<0.5)
- class EasyFEA.Models.Elastic.Orthotropic(
- dim: int,
- E1: float,
- E2: float,
- E3: float,
- G23: float,
- G13: float,
- G12: float,
- v23: float,
- v13: float,
- v12: float,
- axis_1: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (1, 0, 0),
- axis_2: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 1, 0),
- planeStress: bool = True,
- thickness: float = 1.0,
Bases:
_ElasticOrthotropic Linearized Elastic material.
- Walpole_Decomposition() tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- _Behavior(dim: int | None = None)[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices in Kelvin Mandel notation.
In 2D:#
C -> C : Epsilon = Sigma [Sxx Syy sqrt(2)*Sxy]
S -> S : Sigma = Epsilon [Exx Eyy sqrt(2)*Exy]
In 3D:#
C -> C : Epsilon = Sigma [Sxx Syy Szz sqrt(2)*Syz sqrt(2)*Sxz sqrt(2)*Sxy]
S -> S : Sigma = Epsilon [Exx Eyy Ezz sqrt(2)*Eyz sqrt(2)*Exz sqrt(2)*Exy]
- _Update() None[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices. in Kelvin Mandel notation
- E1: float#
Young’s modulus along axis_1.
- E2: float#
Young’s modulus along axis_2.
- E3: float#
Young’s modulus along axis_3.
- G12: float#
Shear modulus in the 1-2 plane.
- G13: float#
Shear modulus in the 1-3 plane.
- G23: float#
Shear modulus in the 2-3 plane.
- v12: float#
Poisson’s ratio for transverse strain along the axis_2 when stressed along the axis_1.
- v13: float#
Poisson’s ratio for transverse strain along the axis_3 when stressed along the axis_1.
- v23: float#
Poisson’s ratio for transverse strain along the axis_3 when stressed along the axis_2.
- class EasyFEA.Models.Elastic.TransverselyIsotropic(
- dim: int,
- El: float,
- Et: float,
- Gl: float,
- vl: float,
- vt: float,
- axis_l: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (1, 0, 0),
- axis_t: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 1, 0),
- planeStress: bool = True,
- thickness: float = 1.0,
Bases:
_ElasticTransversely Isotropic Linearized Elastic material.
- Walpole_Decomposition() tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- _Behavior(dim: int | None = None)[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices in Kelvin Mandel notation.
In 2D:#
C -> C : Epsilon = Sigma [Sxx Syy sqrt(2)*Sxy]
S -> S : Sigma = Epsilon [Exx Eyy sqrt(2)*Exy]
In 3D:#
C -> C : Epsilon = Sigma [Sxx Syy Szz sqrt(2)*Syz sqrt(2)*Sxz sqrt(2)*Sxy]
S -> S : Sigma = Epsilon [Exx Eyy Ezz sqrt(2)*Eyz sqrt(2)*Exz sqrt(2)*Exy]
- _Update() None[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices. in Kelvin Mandel notation
- El: float#
Longitudinal Young’s modulus.
- Et: float#
Transverse Young’s modulus.
- Gl: float#
Longitudinal shear modulus.
- vl: float#
Longitudinal Poisson’s ratio (-1<vl<0.5).
- vt: float#
Transverse Poisson ratio (-1<vt<1)
- class EasyFEA.Models.Elastic._Elastic(dim: int, thickness: float, planeStress: bool)[source]#
Bases:
_IModel,ABCLinearized Elasticity material.
ElasIsot, ElasIsotTrans and ElasAnisot inherit from _Elas class.
- Get_sqrt_C_S() tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#
Returns the Matrix square root of C and S.
- abstractmethod Walpole_Decomposition() tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- _Apply_basis_transformation(
- dim: int,
- material_cM: ndarray[tuple[Any, ...], dtype[floating]],
- material_sM: ndarray[tuple[Any, ...], dtype[floating]],
- axis_1: ndarray[tuple[Any, ...], dtype[floating]],
- axis_2: ndarray[tuple[Any, ...], dtype[floating]],
Performs a basis transformation from the material’s (1,2,3) coordinate system to the (x,y,z) coordinate system to orient the material in space.
- Parameters:
dim (int) – dimension
material_cM (_types.FloatArray) – stiffness matrix
material_sM (_types.FloatArray) – compliance matrix
axis_1 (_types.FloatArray) – Axis 1
axis_2 (_types.FloatArray) – Axis 2
- Returns:
global_cM, global_sM
- Return type:
tuple[_types.FloatArray, _types.FloatArray]
- abstractmethod _Update() None[source]#
Updates the constitutives laws by updating the C stiffness and S compliance matrices. in Kelvin Mandel notation
- property C: ndarray[tuple[Any, ...], dtype[floating]][source]#
Stifness matrix in Kelvin Mandel notation such that:
In 2D: C -> C: Epsilon = Sigma [Sxx, Syy, sqrt(2)*Sxy]
In 3D: C -> C: Epsilon = Sigma [Sxx, Syy, Szz, sqrt(2)*Syz, sqrt(2)*Sxz, sqrt(2)*Sxy].
(Lame’s law)
- property S: ndarray[tuple[Any, ...], dtype[floating]][source]#
Compliance matrix in Kelvin Mandel notation such that:
In 2D: S -> S : Sigma = Epsilon [Exx, Eyy, sqrt(2)*Exy]
In 3D: S -> S: Sigma = Epsilon [Exx, Eyy, Ezz, sqrt(2)*Eyz, sqrt(2)*Exz, sqrt(2)*Exy].
(Hooke’s law)
- dim: int#
- planeStress: bool#
the model uses plane stress simplification
- thickness: float#
Hyper elastic module.
- class EasyFEA.Models.HyperElastic.CiarletGeymonat(
- dim: int,
- K1: float | ndarray[tuple[Any, ...], dtype[floating]],
- K2: float | ndarray[tuple[Any, ...], dtype[floating]],
- K: float | ndarray[tuple[Any, ...], dtype[floating]] = 0.0,
- thickness=1.0,
Bases:
_HyperElastic- Compute_W(
- hyperElasticState: HyperElasticState,
Computes the quadratic energy W(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
We_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2Wde(
- hyperElasticState: HyperElasticState,
Computes dΣde.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
dΣde_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dWde(
- hyperElasticState: HyperElasticState,
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
FeArray – Σ_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
Σxx, Σyy, Σzz, sqrt(2) Σyz, sqrt(2) Σxz, sqrt(2) Σxy
- K: float#
Bulk modulus
- K1: float#
Kappa1
- K2: float#
Kappa2
- class EasyFEA.Models.HyperElastic.HolzapfelOgden(
- dim: int,
- C0: float,
- C1: float,
- C2: float,
- C3: float,
- C4: float,
- C5: float,
- C6: float,
- C7: float,
- K: float,
- Mu1: float,
- Mu2: float,
- T1: ndarray[tuple[Any, ...], dtype[floating]],
- T2: ndarray[tuple[Any, ...], dtype[floating]],
- ks: float = 100,
- thickness=1.0,
Bases:
_HyperElastic- Compute_W(
- hyperElasticState: HyperElasticState,
Computes the quadratic energy W(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
We_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2Wde(
- hyperElasticState: HyperElasticState,
Computes dΣde.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
dΣde_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dWde(
- hyperElasticState: HyperElasticState,
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
FeArray – Σ_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
Σxx, Σyy, Σzz, sqrt(2) Σyz, sqrt(2) Σxz, sqrt(2) Σxy
- C0: float#
- C1: float#
- C2: float#
- C3: float#
- C4: float#
- C5: float#
- C6: float#
- C7: float#
- K: float#
Bulk modulus
- Mu1: float#
- Mu2: float#
- T1#
direction(s) 1, used for the invariants I4 and I8
- T2#
direction(s) 2, used for the invariants I6 and I8
- class EasyFEA.Models.HyperElastic.HyperElasticState(
- mesh: Mesh,
- u: ndarray[tuple[Any, ...], dtype[floating]],
- matrixType: MatrixType,
Bases:
objectHyperelastic state.
- static _CheckFormat(
- mesh: Mesh,
- u: ndarray[tuple[Any, ...], dtype[floating]],
- matrixType: MatrixType,
- Compute_C() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes the right Cauchy-Green tensor C(u) = F(u)’.F(u)
- Returns:
FeArray – C_e_pg of shape (Ne, pg, 3, 3)
dim = 1
——-
cxx 0 0
0 1 0
0 0 1
dim = 2
——-
cxx cxy 0
cyx cyy 0
0 0 1
dim = 3
——-
cxx cxy cxz
cyx cyy cyz
czx czy czz
- Compute_De() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes De(u) derivative green lagrange matrix.
Returns if dim = 2#
- FeArray.FeArrayALike
D_e_pg of shape (Ne, pg, 3, 4)
[1+dxux, 0, dxuy, 0] # xx
[0, dyux, 0, 1+dyuy] # yy
2**(-1/2) [dyux, 1+dxux, 1+dyuy, dxuy # xy
Returns if dim = 3#
- FeArray.FeArrayALike
D_e_pg of shape (Ne, pg, 6, 9)
[1+dxux, 0, 0, dxuy, 0, 0, dxuz, 0, 0] # xx
[0, dyux, 0, 0, 1+dyuy, 0, 0, dyuz, 0] # yy
[0, 0, dzux, 0, 0, dzuy, 0, 0, 1+dzuz] # zz
2**(-1/2) [0, dzux, dyux, 0, dzuy, 1 + dyuy, 0, 1 + dzuz, dyuz] # yz
2**(-1/2) [dzux, 0, 1 + dxux, dzuy, 0, dxuy, 1 + dzuz, 0, dxuz] # xz
2**(-1/2) [dyux, 1+dxux, 0, 1+dyuy, dxuy, 0, dyuz, dxuz, 0] # xy
- Compute_Epsilon() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes the linearized deformation Epsilon = 1/2 (grad(u)’ + grad(u))
Returns if dim = 2#
- FeArray
Eps_e_pg of shape (Ne, pg, 3)
[xx, yy, 2**(-1/2) xy]
Returns if dim = 3#
- FeArray
Eps_e_pg of shape (Ne, pg, 6)
[xx, yy, zz, 2**(-1/2) yz, 2**(-1/2) xz, 2**(-1/2) xy]
- Compute_F() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes the deformation gradient F(u) = I + grad(u)
- Returns:
FeArray – F(u) of shape (Ne, pg, 3, 3)
dim = 1
——-
1+dxux 0 0
0 1 0
0 0 1
dim = 2
——-
1+dxux dyux 0
dxuy 1+dyuy 0
0 0 1
dim = 3
——-
1+dxux dyux dzux
dxuy 1+dyuy dzuy
dxuz dyuz 1+dzuz
- Compute_GreenLagrange() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes the Green-Lagrange deformation E = 1/2 (C - I)
- Returns:
E_e_pg of shape (Ne, pg, dim, dim)
- Return type:
- Compute_I1() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes I1(u)
- Returns:
I1_e_pg of shape (Ne, pg)
- Return type:
- Compute_I2() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes I2(u)
- Returns:
I2_e_pg of shape (Ne, pg)
- Return type:
- Compute_I3() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes I3(u)
- Returns:
I3_e_pg of shape (Ne, pg)
- Return type:
- Compute_I4(
- T: ndarray[tuple[Any, ...], dtype[floating]],
Computes I4(u)
- Parameters:
T (_types.FloatArray) – direction(s)
- Returns:
I4_e_pg of shape (Ne, pg)
- Return type:
- Compute_I6(
- T: ndarray[tuple[Any, ...], dtype[floating]],
Computes I6(u)
- Parameters:
T (_types.FloatArray) – direction(s)
- Returns:
I6_e_pg of shape (Ne, pg)
- Return type:
- Compute_I8(
- T1: ndarray[tuple[Any, ...], dtype[floating]],
- T2: ndarray[tuple[Any, ...], dtype[floating]],
Computes I8(u)
- Parameters:
T1 (_types.FloatArray) – direction(s) 1
T2 (_types.FloatArray) – direction(s) 2
- Returns:
I8_e_pg of shape (Ne, pg)
- Return type:
- Compute_J() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes the deformation gradient J = det(F)
- Returns:
J_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2I1dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes d2I1dC(u)
- Returns:
d2I1dC of shape (d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_d2I2dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes d2I2dC(u)
- Returns:
d2I2dC of shape (d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_d2I3dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes d2I3dC(u)
- Returns:
d2I3dC_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_d2I4dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes d2I4dC(u)
- Returns:
d2I4dC of shape (d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_d2I6dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes d2I6dC(u)
- Returns:
d2I6dC of shape (d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_d2I8dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes d2I8dC(u)
- Returns:
d2I8dC of shape (d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dI1dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes dI1dC(u)
- Returns:
dI1dC of shape (d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dI2dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes dI2dC(u)
- Returns:
dI2dC_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dI3dC() FeArray | ndarray[tuple[Any, ...], dtype[Any]][source]#
Computes dI3dC(u)
- Returns:
dI3dC_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dI4dC(
- T: ndarray[tuple[Any, ...], dtype[floating]],
Computes dI4dC(u)
- Parameters:
T (_types.FloatArray) – direction(s)
- Returns:
dI4dC_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dI6dC(
- T: ndarray[tuple[Any, ...], dtype[floating]],
Computes dI6dC(u)
- Parameters:
T (_types.FloatArray) – direction(s)
- Returns:
dI6dC_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dI8dC(
- T1: ndarray[tuple[Any, ...], dtype[floating]],
- T2: ndarray[tuple[Any, ...], dtype[floating]],
Computes dI8dC(u)
- Parameters:
T1 (_types.FloatArray) – direction(s) 1
T2 (_types.FloatArray) – direction(s) 2
- Returns:
dI8dC_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- _Compute_Anisotropic_Invariants(
- T1: ndarray[tuple[Any, ...], dtype[floating]],
- T2: ndarray[tuple[Any, ...], dtype[floating]],
- _Compute_Anisotropic_Invariants_First_Derivatives(
- T1: ndarray[tuple[Any, ...], dtype[floating]],
- T2: ndarray[tuple[Any, ...], dtype[floating]],
- _Compute_C() list[FeArray | ndarray[tuple[Any, ...], dtype[Any]]][source]#
Computes the right Cauchy-Green tensor components C(u) = F(u)’.F(u)
returns [cxx, cxy, cxz, cyx, cyy, cyz, czx, czy, czz]
- class EasyFEA.Models.HyperElastic.MooneyRivlin(
- dim: int,
- K1: float | ndarray[tuple[Any, ...], dtype[floating]],
- K2: float | ndarray[tuple[Any, ...], dtype[floating]],
- K: float | ndarray[tuple[Any, ...], dtype[floating]] = 0.0,
- thickness=1.0,
Bases:
_HyperElastic- Compute_W(
- hyperElasticState: HyperElasticState,
Computes the quadratic energy W(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
We_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2Wde(
- hyperElasticState: HyperElasticState,
Computes dΣde.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
dΣde_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dWde(
- hyperElasticState: HyperElasticState,
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
FeArray – Σ_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
Σxx, Σyy, Σzz, sqrt(2) Σyz, sqrt(2) Σxz, sqrt(2) Σxy
- K: float#
Bulk modulus
- K1: float#
Kappa1
- K2: float#
Kappa2
- class EasyFEA.Models.HyperElastic.NeoHookean(
- dim: int,
- K: float | ndarray[tuple[Any, ...], dtype[floating]],
- thickness=1.0,
Bases:
_HyperElastic- Compute_W(
- hyperElasticState: HyperElasticState,
Computes the quadratic energy W(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
We_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2Wde(
- hyperElasticState: HyperElasticState,
Computes dΣde.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
dΣde_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dWde(
- hyperElasticState: HyperElasticState,
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
FeArray – Σ_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
Σxx, Σyy, Σzz, sqrt(2) Σyz, sqrt(2) Σxz, sqrt(2) Σxy
- K: float#
Bulk modulus
- class EasyFEA.Models.HyperElastic.SaintVenantKirchhoff(
- dim: int,
- lmbda: float | ndarray[tuple[Any, ...], dtype[floating]],
- mu: float | ndarray[tuple[Any, ...], dtype[floating]],
- K: float | ndarray[tuple[Any, ...], dtype[floating]] = 0.0,
- thickness=1.0,
Bases:
_HyperElastic- Compute_W(
- hyperElasticState: HyperElasticState,
Computes the quadratic energy W(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
We_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2Wde(
- hyperElasticState: HyperElasticState,
Computes dΣde.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
dΣde_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- Compute_dWde(
- hyperElasticState: HyperElasticState,
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
FeArray – Σ_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
Σxx, Σyy, Σzz, sqrt(2) Σyz, sqrt(2) Σxz, sqrt(2) Σxy
- K: float#
Bulk modulus
- lmbda: float#
Lame’s first parameter
- mu: float#
Shear modulus
- class EasyFEA.Models.HyperElastic._HyperElastic(dim: int, thickness: float)[source]#
Bases:
_IModel,ABCHyperElastic material.
NeoHookean, MooneyRivlin, SaintVenantKirchhoff and HolzapfelOgden inherit from _HyperElas class.
- Compute_Tangent_and_Residual(
- hyperElasticState: HyperElasticState,
Computes the tangent matrix and the residual vector.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
tangent_e_pg (Ne, pg, d, d), with d = dof_n * mesh.nPe
- Return type:
- abstractmethod Compute_W(
- hyperElasticState: HyperElasticState,
Computes the quadratic energy W(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
We_e_pg of shape (Ne, pg)
- Return type:
- abstractmethod Compute_d2Wde(
- hyperElasticState: HyperElasticState,
Computes dΣde.
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
dΣde_e_pg of shape (Ne, pg, d, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
- Return type:
- abstractmethod Compute_dWde(
- hyperElasticState: HyperElasticState,
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Returns:
FeArray – Σ_e_pg of shape (Ne, pg, d), where d = 1, 3, 6 depending on whether the solution dimension is 1D, 2D, or 3D.
Σxx, Σyy, Σzz, sqrt(2) Σyz, sqrt(2) Σxz, sqrt(2) Σxy
- dim: float#
- thickness: float#
Beam module.
- class EasyFEA.Models.Beam.BeamStructure(beams: list[_Beam])[source]#
Bases:
_IModelBeam structure class.
- Calc_D_e_pg(
- groupElem: _GroupElem,
- matrixType: MatrixType = MatrixType.beam,
Returns a matrix characterizing the beams’s stiffness behavior.
- Calc_M_e_pg(
- groupElem: _GroupElem,
Returns a matrix characterizing the beams’s stiffness behavior.
- Get_axis_e(
- groupElem: _GroupElem,
Returns the fiber and cross bar vertical axis on every elements.
return xAxis_e, yAxis_e
- property dim: int[source]#
model dimensions
1D -> tension compression
2D -> tension compression + bending + flexion
3D -> all
- class EasyFEA.Models.Beam.Isotropic(
- dim: int,
- line: Line,
- section: Mesh,
- E: float,
- v: float,
- yAxis: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 1, 0),
Bases:
_BeamIsotropic elastic beam.
- Get_D(
- useTimoshenko: bool = False,
Returns a matrix characterizing the beam’s stiffness behavior.
- Get_M() ndarray[tuple[Any, ...], dtype[floating]][source]#
Returns a matrix characterizing the beam’s mass behavior.
- E: float#
Young’s modulus
- v: float#
Poisson’s ratio
- class EasyFEA.Models.Beam._Beam(
- dim: int,
- line: Line,
- section: Mesh,
- yAxis: ndarray[tuple[Any, ...], dtype[Any]] | Collection[int | float] = (0, 1, 0),
Bases:
_IModelBeam class model.
- abstractmethod Get_D(
- useTimoshenko: bool = False,
Returns a matrix characterizing the beam’s stiffness behavior.
- abstractmethod Get_M() ndarray[tuple[Any, ...], dtype[floating]][source]#
Returns a matrix characterizing the beam’s mass behavior.
- _Calc_P() ndarray[tuple[Any, ...], dtype[floating]][source]#
P matrix use to transform beam coordinates to global coordinates.
[ix, jx, kx
iy, jy, ky
iz, jz, kz]
coord(x,y,z) = P • coord(i,j,k)
- _Get_shear_correction_factor(axis: str = 'y') float[source]#
Cowper’s (1966) shear correction factor k for the cross-section, evaluated at Poisson’s ratio ν = 0.
Used to populate
self._ky/self._kzin__init__— those attributes are whatIsotropic.Get_Dactually reads, so users can override them with any value (Pure Jouravski, a textbook table value, Cowper-with-ν=v, …) without touching this helper.Computes k by solving one Saint-Venant Poisson problem on the 2-D section mesh:
∇²φ = -s in the section S ∂φ/∂n = 0 on the boundary ∂S
where s is the section coordinate along the shear direction:
axis="y"→ s = y (gives k_y, pairs with Iz)axis="z"→ s = x (gives k_z, pairs with Iy). The rigid-body (constant) mode of the Neumann problem is fixed by clamping one node to 0; k below is translation-invariant.- Using the Neumann BC, the energy identity gives
uᵀ · f = ∫_S |∇φ|² dS = ∫_S s · φ dS
- and Cowper at ν = 0 reduces to
k = I² / (A · uᵀ · f)
with I = Iz for axis=”y”, I = Iy for axis=”z”.
- Reference values (any mesh size, any aspect ratio):
rectangle → 5/6 ≈ 0.8333 circle → 6/7 ≈ 0.8571 I-beam, hollow tube, channel, … → whatever the geometry yields
- Note vs. the textbook “Pure Jouravski” value:
The two formulas agree on the rectangle (both 5/6).
For curved boundaries they differ: Pure Jouravski’s 1-D τ_yz = V·Q/(I·b), τ_xz = 0 assumption gives 9/10 for a circle; this Poisson PDE returns 6/7 because it accounts for the non-zero τ_xz that the curved boundary forces. 6/7 is the physically correct ν = 0 value (matches 3-D elasticity).
Cowper-with-ν ≠ 0 needs a different PDE (non-zero Neumann boundary term involving ν); not implemented here. If you need a specific k value other than Cowper-ν=0 — Pure Jouravski, Cowper-with-ν, a steel-table value — set
beam._ky/beam._kzdirectly.
- _ky#
Shear correction factor for transverse shear in the beam’s y direction (paired with Iz bending in
Isotropic.Get_D).Populated in
_Beam.__init__with the Cowper-ν=0 value returned by_Get_shear_correction_factor("y"). Override directly to use a different k:# Pure Jouravski (textbook): same as default for rectangle, 9/10 for circle beam._ky = 9 / 10
# Cowper-with-ν for a rectangular section (k = 10(1+ν)/(12+11ν)) beam._ky = 10 * (1 + beam.v) / (12 + 11 * beam.v)
# Cowper-with-ν for a circular section (k = 6(1+ν)/(7+6ν)) beam._ky = 6 * (1 + beam.v) / (7 + 6 * beam.v)
# Value from a steel-section table, etc. beam._ky = 0.49
- _kz#
Shear correction factor for transverse shear in the beam’s z direction (paired with Iy bending in
Isotropic.Get_D). 3-D beams only — unused in 1-D / 2-D simulations.Populated in
_Beam.__init__with the Cowper-ν=0 value returned by_Get_shear_correction_factor("z"). Override directly to use a different k — same Cowper-with-ν / Pure Jouravski / table-value patterns as_ky, just applied to the z direction (so e.g. for a rectangle bxh, the paired bending inertia is Iy = h·b³/12 and the formula factor refers to the b direction).
- property dof_n: int[source]#
Degrees of freedom per node 1D -> [u1, … , un]
2D -> [u1, v1, rz1, …, un, vn, rzn]
3D -> [u1, v1, w1, rx1, ry1, rz1, …, u2, v2, w2, rx2, ry2, rz2]