simulations#

The EasyFEA/simulations/ module in EasyFEA provides essential tools for creating and managing simulations. These simulations are built using a Mesh and a _IModel (material).

With this module, you can construct:

How to create new simulations in EasyFEA ?#

To create new simulation classes, you can take inspiration from existing implementations. Make sure to follow the _Simu interface. The ThermalSimu class is relatively simple and can serve as a good starting point. See EasyFEA/simulations/_thermal.py for more details.

Detailed simulations API#

class EasyFEA.simulations.AlgoType(value)[source]#

Bases: str, Enum

static Get_Hyperbolic_Types()[source]#
Return type:

list[str]

static Get_Hyperbolic_and_Parabolic_Types()[source]#
Return type:

list[str]

elliptic = 'elliptic'#

Solve K u = F

hht = 'hht'#

Solve K u_np1ma + C v_np1ma + M a_np1ma = F_np1ma

midpoint = 'midpoint'#

Solve K u_np1/2 + C v_np1/2 + M a_np1/2 = F_np1/2

newmark = 'newmark'#

Solve K u_np1 + C v_np1 + M a_np1 = F_np1

u_np1 = u_n + dt v_n + dt^2/2 ((1 - 2 B) a_n + 2 B a_np1)

v_np1 = v_n + dt ((1 - gamma) a_n + gamma a_np1)

parabolic = 'parabolic'#

Solve K u_npa + C v_npa = F_npa

EasyFEA.simulations.Load_Simu(folder, filename='simulation')[source]#

Loads the simulation from the specified folder.

Parameters:
  • folder (str) – simulation’s folder.

  • filename (str, optional) – The simualtion’s name, by default “simulation”.

Returns:

The loaded simulation.

Return type:

_Simu

EasyFEA.simulations.Load_pickle(folder, filename)[source]#

Returns folder/filename.pickle object.

class EasyFEA.simulations.ResolType(value)[source]#

Bases: str, Enum

r1 = '1'#

xi = inv(Aii) * (bi - Aic * xc)

r2 = '2'#

Lagrange multipliers

r3 = '3'#

Penality

EasyFEA.simulations.Save_pickle(obj, folder, filename)[source]#

Saves the object in folder/filename.pickle.

Return type:

None

EasyFEA.simulations.Solve_simu(simu, problemType)[source]#

Solving the simulation’s problem according to the resolution type.

class EasyFEA.simulations._Simu(mesh, model, verbosity=True, useIterativeSolvers=True)[source]#

Bases: _IObserver, Updatable, ABC

The following classes inherit from the parent class _Simu:
  • ElasticSimu

  • HyperElasticSimu

  • PhaseFieldSimu

  • BeamSimu

  • ThermalSimu

  • WeakFormSimu

To create new simulation classes, you can take inspiration from existing implementations.

Make sure to follow this interface.

The ThermalSimuclass is relatively simple and can serve as a good starting point.

See simulations/_thermal.py for more details.

To use the interface/inheritance, 13 methods need to be defined.

General:#

  • def Get_problemTypes(self) -> list[ModelType]:

  • def Get_unknowns(self, problemType=None) -> list[str]:

  • def Get_dof_n(self, problemType=None) -> int:

These functions provides access to the available unknowns and degrees of freedom (dofs).

Solve:#

  • def Get_x0(self, problemType=None):

  • def Construct_local_matrix_system(self, problemType):

These functions assemble the matrix system K u + C v + M a = F.

Iterations:#

  • def Save_Iter(self) -> None:

  • def Set_Iter(self, index=-1) -> None:

These functions are used to save or load iterations.

Results:#

  • def Results_Available(self) -> list[str]:

  • def Result(self, result: str, nodeValues=True, iter=None) -> float | _types.FloatArray:

  • def Results_Iter_Summary(self) -> tuple[list[int], list[tuple[str, _types.FloatArray]]]:

  • def Results_dict_Energy(self) -> dict[str, float]:

  • def Results_displacement_matrix(self) -> _types.FloatArray:

  • def Results_nodesField_elementsField(self, details=False) -> tuple[list[str], list[str]]:

These functions are used to process the results.

Assembly(problemType)[source]#

Assembles the matrix system K u + C v + M a = F for the given problemType.

Return type:

tuple[csr_matrix, csr_matrix, csr_matrix, csr_matrix]

property Bc_Dirichlet: list[BoundaryCondition]#

Returns a copy of the Dirichlet conditions.

property Bc_Display: list[BoundaryCondition | LagrangeCondition]#

Returns a copy of the boundary conditions for display.

Bc_Init()[source]#

Initializes Dirichlet, Neumann and Lagrange boundary conditions

Return type:

None

property Bc_Lagrange: list[LagrangeCondition]#

Returns a copy of the Lagrange conditions.

property Bc_Neuman: list[BoundaryCondition]#

Returns a copy of the Neumann conditions.

Bc_dofs_Dirichlet(problemType=None)[source]#

Returns dofs related to Dirichlet conditions.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Bc_dofs_Neumann(problemType=None)[source]#

Returns dofs related to Neumann conditions.

Return type:

ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]

Bc_dofs_known_unknown(problemType)[source]#

Returns known and unknown dofs.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]], ndarray[tuple[Any, ...], dtype[TypeVar(IntType, bound= Union[integer, int])]]]

Bc_dofs_nodes(nodes, unknowns, problemType=None)[source]#

Returns degrees of freedom associated with the nodes, based on the problem type and unknowns.

Parameters:
  • nodes (_types.IntArray) – nodes.

  • unknowns (list) – unknowns (e.g [“x”,”y”,”rz”])

  • problemType (str) – Problem type.

Returns:

Degrees of freedom.

Return type:

_types.IntArray

Bc_values_Dirichlet(problemType=None)[source]#

Returns dofs values related to Dirichlet conditions.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Bc_values_Neumann(problemType=None)[source]#

Returns dofs values related to Neumann conditions.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Bc_vector_Dirichlet(problemType=None)[source]#

Returns a vector filled with Dirichlet boundary conditions values.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Bc_vector_Neumann(problemType=None)[source]#

Returns a vector filled with Neuman boundary conditions values.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

abstractmethod Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Return type:

tuple[Optional[FeArray], Optional[FeArray], Optional[FeArray], Optional[FeArray]]

Get_K_C_M_F(problemType=None)[source]#

Returns the assembled matrices of K u + C v + M a = F for the given problem.

Return type:

tuple[csr_matrix, csr_matrix, csr_matrix, csr_matrix]

Get_contact(masterMesh, slaveNodes=None, masterNodes=None)[source]#

Returns the simulation nodes detected in the master mesh with the associated displacement matrix to the interface.

Parameters:
  • masterMesh (Mesh) – master mesh

  • slavenodes (_types.IntArray, optional) – slave nodes, by default None

  • masternodes (_types.IntArray, optional) – master nodes, by default None

Returns:

nodes, displacementMatrix

Return type:

tuple[_types.IntArray, _types.FloatArray]

abstractmethod Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_lb_ub(problemType)[source]#

Returns the lower bound and upper bound.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]

abstractmethod Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Return type:

list[ModelType]

abstractmethod Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

abstractmethod Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Need_Update(value=True)[source]#

Sets whether the simulation needs to reconstruct matrices K, C, M and F.

property Niter: int#

Number of iterations

abstractmethod Result(option, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float]

abstractmethod Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Get_Bc_Summary()[source]#

Returns the simulation loading summary.

Return type:

str

Results_Get_Iteration_Summary()[source]#

Returns the iteration’s summary.

Return type:

str

abstractmethod Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_Reshape_values(values, nodeValues)[source]#

Reshapes input values based on whether they are stored at nodes or elements.

Parameters:
  • values (_types.FloatArray) – Input values to reshape.

  • nodeValues (bool) – If True, the output will represent values at nodes; if False, values on elements will be derived.

Returns:

Reshaped values on nodes or elements.

Return type:

_types.FloatArray

Raises:

Exception – Raised if unexpected conditions occur during the calculation.

Results_Set_Bc_Summary()[source]#

Sets the simulation loading summary.

Return type:

None

Results_Set_Iteration_Summary()[source]#

Sets the iteration’s summary.

Return type:

None

abstractmethod Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Return type:

dict[str, float]

abstractmethod Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

abstractmethod Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Return type:

tuple[list[str], list[str]]

Save(folder, filename='simulation', additionalInfos='')[source]#

Saves the simulation and its summary in the folder. Saves the simulation as ‘filename.pickle’.

Return type:

None

abstractmethod Save_Iter()[source]#

Saves iteration results in _results.

Return type:

dict[str, Any]

abstractmethod Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

Return type:

dict

Solve()[source]#

Computes the solution field for the current boundary conditions.

Returns:

The solution of the simulation.

Return type:

_types.FloatArray

Solver_Set_Elliptic_Algorithm()[source]#

Sets the algorithm’s resolution properties for an elliptic problem.

Used to solve K u = F.

Return type:

None

Solver_Set_Hyperbolic_Algorithm(dt, beta=0.25, gamma=0.5, algo=AlgoType.newmark, alpha=0.5)[source]#

Sets the algorithm’s resolution properties for a Hyperbolic problem.

Used to solve K u + C v + M a = F.

Parameters:
  • dt (float) – The time increment.

  • beta (float, optional) – The coefficient beta, by default 1/4.

  • gamma (float, optional) – The coefficient gamma, by default 1/2.

  • algo (AlgoType, optional) –

    Algo used to solve K u + C v + M a = F, by default AlgoType.newmark

    K u_np1 + C v_np1 + M a_np1 = F_np1.

  • alpha (float, optional) – The coefficient alpha, by default 1/2.

Return type:

None

Solver_Set_Parabolic_Algorithm(dt, alpha=0.5)[source]#

Sets the algorithm’s resolution properties for a parabolic problem.

Used to solve K u_np1 + C v_np1 = F_np1 with:

u_np1 = u_n + dt v_npa

Parameters:
  • dt (float) – The time increment.

  • alpha (float, optional) –

    The alpha criterion, by default 1/2

    • 0 -> Forward Euler

    • 1 -> Backward Euler

    • 1/2 -> midpoint

Return type:

None

add_dirichlet(nodes, values, unknowns, problemType=None, description='')[source]#

Adds Dirichlet’s boundary conditions.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contains floats, arrays or functions or functions.

    e.g [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    The functions use the x, y and z nodes coordinates.

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

add_lineLoad(nodes, values, unknowns, problemType=None, description='')[source]#

Adds a linear load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

add_neumann(nodes, values, unknowns, problemType=None, description='')[source]#

Adds Neumann’s boundary conditions.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contains floats, arrays or functions or functions.

    e.g [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    The functions use the x, y and z nodes coordinates.

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

add_pressureLoad(nodes, magnitude, problemType=None, description='')[source]#

Adds a pressure.

Parameters:
  • nodes (_types.IntArray) – nodes. (must belong to the edge of the mesh.)

  • magnitude (float) – pressure magnitude

  • problemType (str, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

add_surfLoad(nodes, values, unknowns, problemType=None, description='')[source]#

Adds a surface load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

add_volumeLoad(nodes, values, unknowns, problemType=None, description='')[source]#

Adds a volumetric load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

property algo: AlgoType#

The algorithm used to solve the problem.

(elliptic, parabolic, hyperbolic) see:

  • Solver_Set_Elliptic_Algorithm()

K u = F - Solver_Set_Parabolic_Algorithm()

K u + C v = F - Solver_Set_Hyperbolic_Algorithm()

K u + C v + M a = F

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

Center of mass / barycenter / inertia center

property dim: int#

simulation’s dimension

property isNonLinear#

Returns whether the simulation is non linear.

property mass: float#
property mesh: Mesh#

simulation’s mesh.

property model: _IModel#

model used

property problemType: ModelType#

Get the simulation problem type.

property results: list[dict]#

Returns a copy of the list of dictionary containing the results from each iteration.

rho#

mass density

property solver: str#

Solver used to solve the simulation.

property useIterativeSolvers: bool#

Iterative solvers can be used.

This module contains all available simulation classes.

class EasyFEA.Simulations.BeamSimu(mesh, model, verbosity=False, useIterativeSolvers=True)[source]#
Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Return type:

list[ModelType]

Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Result(result, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float]

Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Get_Iteration_Summary()[source]#

Returns the iteration’s summary.

Return type:

str

Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Return type:

dict[str, float]

Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Return type:

tuple[list[str], list[str]]

Save_Iter()[source]#

Saves iteration results in _results.

Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

Return type:

dict

add_connection(nodes, unknowns, description)[source]#

Connects beams together in the specified unknowns.

Parameters:
  • nodes (_types.IntArray) – nodes

  • unknowns (list[str]) – unknowns

  • description (str) – description

add_connection_fixed(nodes, description='Fixed')[source]#

Adds a fixed connection.

Parameters:
  • nodes (_types.IntArray) – nodes

  • description (str, optional) – description, by default “Fixed”

add_connection_hinged(nodes, unknowns=[''], description='Hinged')[source]#

Adds a hinged connection.

Parameters:
  • nodes (_types.IntArray) – nodes

  • unknowns (list, optional) – unknowns, by default [‘’]

  • description (str, optional) – description, by default “Hinged”

add_surfLoad(nodes, values, unknowns, problemType=None, description='')[source]#

Adds a surface load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

add_volumeLoad(nodes, values, unknowns, problemType=None, description='')[source]#

Adds a volumetric load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

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

Center of mass / barycenter / inertia center

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

Displacement vector field.

1D [uxi, …]

2D [uxi, uyi, rzi, …]

3D [uxi, uyi, uzi, rxi, ryi, rzi, …]

property structure: BeamStructure#

Beam structure.

class EasyFEA.Simulations.ElasticSimu(mesh, model, verbosity=False, useIterativeSolvers=True)[source]#
Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Return type:

list[ModelType]

Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Result(result, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float]

Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Get_Iteration_Summary()[source]#

Returns the iteration’s summary.

Return type:

str

Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Return type:

dict[str, float]

Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Return type:

tuple[list[str], list[str]]

Save_Iter()[source]#

Saves iteration results in _results.

Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

Return type:

dict

Set_Rayleigh_Damping_Coefs(coefM=0.0, coefK=0.0)[source]#

Sets damping coefficients ( C = coefK * K + coefM * M ).

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

Acceleration vector field.

2D [axi, ayi, …]

3D [axi, ayi, azi, …]

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

Displacement vector field.

2D [uxi, uyi, …]

3D [uxi, uyi, uzi, …]

property material: _Elas#

elastic material

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

Velocity vector field.

2D [vxi, vyi, …]

3D [vxi, vyi, vzi, …]

class EasyFEA.Simulations.HyperElasticSimu(mesh, model, tolConv=1e-05, maxIter=20, verbosity=False, useIterativeSolvers=True)[source]#
Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Result(result, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float]

Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Save_Iter()[source]#

Saves iteration results in _results.

Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

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

Displacement vector field.

[uxi, uyi, uzi, …]

property material: _HyperElas#

hyperelastic material

class EasyFEA.Simulations.PhaseFieldSimu(mesh, model, verbosity=False, useIterativeSolvers=True)[source]#
Bc_dofs_nodes(nodes, unknowns, problemType=ModelType.elastic)[source]#

Returns degrees of freedom associated with the nodes, based on the problem type and unknowns.

Parameters:
  • nodes (_types.IntArray) – nodes.

  • unknowns (list) – unknowns (e.g [“x”,”y”,”rz”])

  • problemType (str) – Problem type.

Returns:

Degrees of freedom.

Return type:

_types.IntArray

Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Get_K_C_M_F(problemType=None)[source]#

Returns the assembled matrices of K u + C v + M a = F for the given problem.

Return type:

tuple[csr_matrix, csr_matrix, csr_matrix, csr_matrix]

Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_lb_ub(problemType)[source]#

Returns the lower bound and upper bound.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]

Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Return type:

list[ModelType]

Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Need_Update(value=True)[source]#

Sets whether the simulation needs to reconstruct matrices K, C, M and F.

Return type:

None

Result(result, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float, None]

Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Get_Bc_Summary()[source]#

Returns the simulation loading summary.

Return type:

str

Results_Get_Iteration_Summary()[source]#

Returns the iteration’s summary.

Return type:

str

Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_Set_Bc_Summary(config)[source]#

Sets the simulation loading summary.

Results_Set_Iteration_Summary(iter, load, unitLoad, percentage=0.0, remove=False)[source]#

Creates the iteration summary for the damage problem.

Parameters:
  • iter (int) – iteration

  • load (float) – loading

  • unitLoad (str) – loading unit

  • percentage (float, optional) – percentage of simualtion performed, by default 0.0

  • remove (bool, optional) – removes line from terminal after display, by default False

Return type:

str

Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Return type:

dict[str, float]

Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Return type:

tuple[list[str], list[str]]

Save_Iter()[source]#

Saves iteration results in _results.

Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

Return type:

dict

Solve(tolConv=1.0, maxIter=500, convOption=2)[source]#

Solves the iterative damage problem using the staggered scheme.

Parameters:
  • tolConv (float, optional) – threshold used to check convergence (𝜖), by default 1.0

  • maxIter (int, optional) – Maximum iterations for convergence, by default 500

  • convOption (int, optional) –

    • 0 -> convergence on damage np.max(np.abs(d_np1-dk)) equivalent normInf(d_np1-dk)

    • 1 -> convergence on crack energy np.abs(psi_crack_n - psi_crack_np1)/psi_crack_np1

    • 2 -> convergence on total energy np.abs(psi_tot_n - psi_tot_np1)/psi_tot_np1

    eq (39) Ambati 2015 10.1007/s00466-014-1109-y

    • 3 -> (convD <= tolConv) and (convU <= tolConv*0.999)

    eq (25) Pech 2022 10.1016/j.engfracmech.2022.108591

Returns:

u_np1, d_np1, Ku, converged

such that:

  • u_np1: displacement vector field

  • d_np1: damage scalar field

  • Ku: displacement stiffness matrix

  • converged: the solution has converged

Return type:

_types.FloatArray, _types.FloatArray, csr_matrix, bool

add_dirichlet(nodes, values, unknowns, problemType=ModelType.elastic, description='')[source]#

Adds Dirichlet’s boundary conditions.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contains floats, arrays or functions or functions.

    e.g [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    The functions use the x, y and z nodes coordinates.

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

add_lineLoad(nodes, values, unknowns, problemType=ModelType.elastic, description='')[source]#

Adds a linear load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

add_neumann(nodes, values, unknowns, problemType=ModelType.elastic, description='')[source]#

Adds Neumann’s boundary conditions.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contains floats, arrays or functions or functions.

    e.g [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    The functions use the x, y and z nodes coordinates.

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

add_pressureLoad(nodes, magnitude, problemType=ModelType.elastic, description='')[source]#

Adds a pressure.

Parameters:
  • nodes (_types.IntArray) – nodes. (must belong to the edge of the mesh.)

  • magnitude (float) – pressure magnitude

  • problemType (str, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

Return type:

None

add_surfLoad(nodes, values, unknowns, problemType=ModelType.elastic, description='')[source]#

Adds a surface load.

Parameters:
  • nodes (_types.IntArray) – nodes

  • values (list) –

    list of values that can contain floats, arrays or functions or lambda functions.

    e.g = [10, lambda x,y,z: 10*x - 20*y + x*z, _types.FloatArray]

    functions use x, y and z integration points coordinates (x,y,z are in this case arrays of dim (e,p))

    Please note that the functions must take 3 input parameters in the order x, y, z, whether the problem is 1D, 2D or 3D.

  • unknowns (list[str]) – unknowns where values will be applied (e.g [‘y’, ‘x’])

  • problemType (ModelType, optional) – problem type, if not specified, we take the basic problem of the problem

  • description (str, optional) – Description of the condition, by default “”.

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

Damage scalar field.

[di, …]

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

Displacement vector field.

2D [uxi, uyi, …]

3D [uxi, uyi, uzi, …]

property needUpdate: bool#

The object needs to be updated.

property phaseFieldModel: PhaseField#

damage model

class EasyFEA.Simulations.ThermalSimu(mesh, model, verbosity=False, useIterativeSolvers=True)[source]#
Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Return type:

list[ModelType]

Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Result(result, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float]

Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Return type:

dict[str, float]

Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Return type:

tuple[list[str], list[str]]

Save_Iter()[source]#

Saves iteration results in _results.

Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

Return type:

dict

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

Scalar temperature field.

[ti, ….]

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

Time derivative of the scalar temperature field.

[d(ti)/dt, ….]

property thermalModel: Thermal#

Thermal simulation model.

class EasyFEA.Simulations.WeakFormSimu(mesh, model, isNonLinear=False, tolConv=1e-05, maxIter=20, verbosity=False, useIterativeSolvers=True)[source]#
Construct_local_matrix_system(problemType)[source]#

Construct the local matrix system K u + C v + M a = F for the given problem.

Get_dof_n(problemType=None)[source]#

Returns the number of degrees of freedom per node.

Return type:

int

Get_problemTypes()[source]#

Returns the problem types available through the simulation.

Return type:

list[ModelType]

Get_unknowns(problemType=None)[source]#

Returns a list of unknowns available in the simulation.

Return type:

list[str]

Get_x0(problemType=None)[source]#

Returns the solution from the previous iteration.

Result(result, nodeValues=True, iter=None)[source]#

Returns the result. Use Results_Available() to know the available results.

Return type:

Union[ndarray[tuple[Any, ...], dtype[floating]], float]

Results_Available()[source]#

Returns a list of available results in the simulation.

Return type:

list[str]

Results_Iter_Summary()[source]#

Returns the values to be displayed in Plot_Iter_Summary.

Return type:

tuple[list[int], list[tuple[str, ndarray[tuple[Any, ...], dtype[floating]]]]]

Results_dict_Energy()[source]#

Returns a dictionary containing the names and values of the calculated energies.

Return type:

dict[str, float]

Results_displacement_matrix()[source]#

Returns displacements as a matrix [dx, dy, dz] (Nn,3).

Return type:

ndarray[tuple[Any, ...], dtype[floating]]

Results_nodeFields_elementFields(details=False)[source]#

Returns lists of nodesFields and elementsFields displayed in paraview.

Return type:

tuple[list[str], list[str]]

Save_Iter()[source]#

Saves iteration results in _results.

Set_Iter(iter=-1, resetAll=False)[source]#

Sets the simulation to the specified iteration (usually the last one) and then reset the required variables if necessary (resetAll).

Returns the simulation results dictionary.

Return type:

dict

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

node field a = d2udt2

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

node field u.

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

node field v = dudt

property weakForms: WeakForms#

Weak form manager.