models#
The EasyFEA/models/ module in EasyFEA provides essential tools for creating and managing models.
These models are used to build _Simu instances and mainly contain material parameters.
With this module, you can construct:
Linear elastic materials, such as
ElasIsot,ElasIsotTrans,ElasOrthotropic, andElasAnisot.Nonlinear hyperelastic materials, such as
NeoHookean,MooneyRivlin,SaintVenantKirchhoff, andHolzapfelOgden.Elastic beams with
BeamElasIsot.Phase-field materials with
PhaseField.Thermal materials with
Thermal.Weark forms with
WeakForms.
Detailed materials API#
Module containing constitutive laws for materials, including elastic, hyperelastic, damage, thermal, and beam materials.
- EasyFEA.models.Apply_Pmat(P, M, toGlobal=True)[source]#
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
- class EasyFEA.models.BeamElasIsot(dim, line, section, E, v, yAxis=(0, 1, 0))[source]#
Bases:
_BeamIsotropic elastic beam.
-
E:
float# Young’s modulus
- Get_D()[source]#
Returns a matrix characterizing the beam’s stiffness behavior.
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
- property mu: float#
shear modulus (G)
-
v:
float# Poisson’s ratio
-
E:
- class EasyFEA.models.BeamStructure(beams)[source]#
Bases:
_IModelBeam structure class.
- Calc_D_e_pg(groupElem)[source]#
Returns a matrix characterizing the beams’s stiffness behavior.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_axis_e(groupElem)[source]#
Returns the fiber and cross bar vertical axis on every elements.
return xAxis_e, yAxis_e
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[floating]]]
- property areas: list[float]#
beams areas
- property beams: list[_Beam]#
- property dim: int#
model dimensions
1D -> tension compression
2D -> tension compression + bending + flexion
3D -> all
- property dof_n: int#
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]
- property nBeam: int#
Number of beams in the structure
- property thickness: float#
The beam structure can have several beams and therefore different sections.
You need to look at the section of the beam you are interested in.
- class EasyFEA.models.ElasAnisot(dim, C, useVoigtNotation, axis1=(1, 0, 0), axis2=(0, 1, 0), thickness=1.0)[source]#
Bases:
_ElasAnisotropic Linearized Elastic material.
- Set_C(C, useVoigtNotation=True, update_S=True)[source]#
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()[source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[floating]]]
- property axis1: ndarray[tuple[Any, ...], dtype[floating]]#
axis1 vector
- property axis2: ndarray[tuple[Any, ...], dtype[floating]]#
axis2 vector
- class EasyFEA.models.ElasIsot(dim, E=210000.0, v=0.3, planeStress=True, thickness=1.0)[source]#
Bases:
_ElasIsotropic Linearized Elastic material.
-
E:
float# Young’s modulus
- Walpole_Decomposition()[source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[floating]]]
-
v:
float# Poisson’s ratio (-1<v<0.5)
-
E:
- class EasyFEA.models.ElasIsotTrans(dim, El, Et, Gl, vl, vt, axis_l=(1, 0, 0), axis_t=(0, 1, 0), planeStress=True, thickness=1.0)[source]#
Bases:
_ElasTransversely Isotropic Linearized Elastic material.
-
El:
float# Longitudinal Young’s modulus.
-
Et:
float# Transverse Young’s modulus.
-
Gl:
float# Longitudinal shear modulus.
- property Gt: float | ndarray[tuple[Any, ...], dtype[floating]]#
Transverse shear modulus.
- Walpole_Decomposition()[source]#
Walpole’s decomposition in Kelvin Mandel notation such that:
C = sum(ci * Ei).
returns ci, Ei
- Return type:
tuple[ndarray[tuple[Any,...],dtype[floating]],ndarray[tuple[Any,...],dtype[floating]]]
- property axis_l: ndarray[tuple[Any, ...], dtype[floating]]#
Longitudinal axis
- property axis_t: ndarray[tuple[Any, ...], dtype[floating]]#
Transversal axis
- property kt: float | ndarray[tuple[Any, ...], dtype[floating]]#
-
vl:
float# Longitudinal Poisson’s ratio (-1<vl<0.5).
-
vt:
float# Transverse Poisson ratio (-1<vt<1)
-
El:
- EasyFEA.models.Get_Pmat(axis_1, axis_2, useMandel=True)[source]#
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”
- class EasyFEA.models.HolzapfelOgden(dim, C0, C1, C2, C3, C4, C5, C6, C7, K, Mu1, Mu2, T1, T2, ks=100, thickness=1.0)[source]#
Bases:
_HyperElas-
C0:
float#
-
C1:
float#
-
C2:
float#
-
C3:
float#
-
C4:
float#
-
C5:
float#
-
C6:
float#
-
C7:
float#
- Compute_W(hyperElasticState)[source]#
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)[source]#
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)[source]#
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Return 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
-
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
-
C0:
- class EasyFEA.models.HyperElasticState(mesh, u, matrixType)[source]#
Bases:
objectHyperelastic state.
- Compute_C()[source]#
Computes the right Cauchy-Green tensor C(u) = F(u)’.F(u)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]- 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()[source]#
Computes De(u) derivative green lagrange matrix.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
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()[source]#
Computes the linearized deformation Epsilon = 1/2 (grad(u)’ + grad(u))
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
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()[source]#
Computes the deformation gradient F(u) = I + grad(u)
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]- 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()[source]#
Computes the Green-Lagrange deformation E = 1/2 (C - I)
- Returns:
E_e_pg of shape (Ne, pg, dim, dim)
- Return type:
- Compute_I4(T)[source]#
Computes I4(u)
- Parameters:
T (_types.FloatArray) – direction(s)
- Returns:
I4_e_pg of shape (Ne, pg)
- Return type:
- Compute_I6(T)[source]#
Computes I6(u)
- Parameters:
T (_types.FloatArray) – direction(s)
- Returns:
I6_e_pg of shape (Ne, pg)
- Return type:
- Compute_I8(T1, T2)[source]#
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()[source]#
Computes the deformation gradient J = det(F)
- Returns:
J_e_pg of shape (Ne, pg)
- Return type:
- Compute_d2I1dC()[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()[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()[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()[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()[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()[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()[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()[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()[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)[source]#
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)[source]#
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, T2)[source]#
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:
- property mesh#
- EasyFEA.models.KelvinMandel_Matrix(dim, M)[source]#
Apply Kelvin Mandel coefficient to constitutive laws.
- Return type:
ndarray[tuple[Any,...],dtype[floating]]
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]]
- 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.MooneyRivlin(dim, K1, K2, K=0.0, thickness=1.0)[source]#
Bases:
_HyperElas- Compute_W(hyperElasticState)[source]#
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)[source]#
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)[source]#
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Return 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.NeoHookean(dim, K, thickness=1.0)[source]#
Bases:
_HyperElas- Compute_W(hyperElasticState)[source]#
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)[source]#
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)[source]#
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Return 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.PhaseField(material, split, regularization, Gc, l0, solver=SolverType.History, A=None)[source]#
Bases:
_IModelPhaseField class.
- A#
matrix characterizing the weak anisotropy in the crack surface density function
- Calc_C(Epsilon_e_pg, verif=False)[source]#
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)[source]#
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)[source]#
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
-
Gc:
float# critical energy release rate
- Get_f_e_pg(PsiP_e_pg)[source]#
Returns source therm
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_g_e_pg(d_n, mesh, matrixType, k_res=1e-12)[source]#
Returns degradation function
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- Get_r_e_pg(PsiP_e_pg)[source]#
Returns reaction therm
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- static Get_solvers()[source]#
Returns available solvers used to manage crack irreversibility
- Return type:
list[SolverType]
- 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'#
- property c_w#
scaling parameter for accurate dissipation of crack energy
- property dim: int#
model dimension
- property isHeterogeneous: bool#
indicates whether the model has heterogeneous parameters
- property k: float | ndarray[tuple[Any, ...], dtype[floating]]#
get diffusion therm
-
l0:
float# half crack width
- property material: _Elas#
elastic material
-
solver:
SolverType# solver used to manage crack irreversibility
- property thickness: float#
thickness used in the model
- EasyFEA.models.Project_Kelvin(A, orderA=None)[source]#
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, coef=np.float64(1.4142135623730951))[source]#
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- EasyFEA.models.Project_vector_to_matrix(vector, coef=np.float64(1.4142135623730951))[source]#
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- EasyFEA.models.Reshape_variable(variable, Ne, nPg)[source]#
Resizes variable to (Ne, nPg, …) shape.
- Return type:
Union[FeArray,ndarray[tuple[Any,...],dtype[Any]]]
- EasyFEA.models.Result_in_Strain_or_Stress_field(field_e_pg, result, coef=np.float64(1.4142135623730951))[source]#
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
- class EasyFEA.models.SaintVenantKirchhoff(dim, lmbda, mu, K=0.0, thickness=1.0)[source]#
Bases:
_HyperElas- Compute_W(hyperElasticState)[source]#
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)[source]#
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)[source]#
Computes the second Piola-Kirchhoff tensor Σ(u).
- Parameters:
hyperElasticState (HyperElasticState) – Hyperelastic state containing the mesh, the discretized field, and the matrix type.
- Return 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.Thermal(k, c=0.0, thickness=1.0)[source]#
Bases:
_IModelThermal class.
- c#
mass heat capacity [J K^-1 kg^-1]
-
dim:
int#
- property isHeterogeneous: bool#
indicates whether the model has heterogeneous parameters
- k#
thermal conductivity [W m^-1]
-
thickness:
float#
- class EasyFEA.models.WeakForms(field, computeK, computeC=None, computeM=None, computeF=None, thickness=1.0)[source]#
Bases:
_IModelClass responsible for computing the finite element matrices used in the system ( K u + C v + M a = F ).
- property computeC: BiLinearForm | None#
Function used to build damping matrix C from ( K u + C v + M a = F ).
- property computeF: LinearForm | None#
Function used to build force vector F from ( K u + C v + M a = F ).
- property computeK: BiLinearForm | None#
Function used to build stiffness matrix K from ( K u + C v + M a = F ).
- property computeM: BiLinearForm | None#
Function used to build mass matrix M from ( K u + C v + M a = F ).
- property dim: int#
model dimension
-
thickness:
float#