FEM package#

FEM implementation for N dimensional and M variables per node.

FEM.importJSON(json_file, aditional, **kargs)#

Subpackages#

Submodules#

FEM.Core module#

Define the structure of a finite element problem. This is the parent class of individual equation classes.

The individual children classes must implement the method for calculating the element matrices and post processing.

class FEM.Core.Core(geometry: Geometry, solver: Lineal | NonLinealSolver | None = None, sparse: bool = False, verbose: bool = False, name='')#

Bases: object

Create the Finite Element problem.

Parameters:
  • geometry (Geometry) – Input geometry. The geometry must contain the elements, and the border conditions.

  • class. (You can create the geometry of the problem using the Geometry) –

  • solver (Union[Lineal, NonLinealSolver], optional) – Finite Element solver. If not provided, Lineal solver is used.

  • sparse (bool, optional) – To use sparse matrix formulation. Defaults to False

  • verbose (bool, optional) – To print console messages and progress bars. Defaults to False.

borderConditions() None#

Assign border conditions to the system. The border conditios are assigned in this order:

  1. Natural border conditions

  2. Essential border conditions

This ensures that in a node with 2 border conditions the essential border conditions will be applied.

condensedSystem() None#

Assign border conditions to the system and modifies the matrices to condensed mode The border conditios are assigned in this order:

  1. Natural border conditions

  2. Essential border conditions

This ensures that in a node with 2 border conditions the essential border conditions will be applied.

description()#

Generates the problem description for loggin porpuses

elementMatrices() None#

Calculate the element matrices

ensembling() None#

Ensembling of equation system. This method use the element gdl and the element matrices. The element matrices degrees of fredom must match the dimension of the element gdl. For m>1 variables per node, the gdl will be flattened. This ensure that the element matrices will always be a 2-D Numpy Array.

exportJSON(filename: str | None = None)#
postProcess() None#

Post process the solution

profile() None#

Create a profile for a 3D or 2D problem.

restartMatrix() None#

Sets all model matrices and vectors to 0 state

solve(plot: bool = True, **kargs) None#

A series of Finite Element steps

Parameters:
  • plot (bool, optional) – To post process the solution. Defaults to True.

  • **kargs – Solver specific parameters.

solveES(**kargs) None#

Solve the finite element problem

solveFromArray(solution: ndarray, plot: bool = True, **kargs) None#

Load a solution array to the problem.

Parameters:
  • solution (np.ndarray) – Solution vertical array with shape (self.ngdl,1)

  • plot (bool, optional) – To post process the solution. Defaults to True.

solveFromFile(file: str, plot: bool = True, **kargs) None#

Load a solution file and show the post process for a given geometry

Parameters:
  • file (str) – Path to the previously generated solution file.

  • plot (bool, optional) – To post process the solution. Defaults to True.

class FEM.Core.CoreHiperbolic(geometry: Geometry, solver: Solver | None = None, verbose: bool = False, name='')#

Bases: CoreTransient

docstring for CoreParabolic

apply_initial_condition()#
set_initial_condition(u0: list[float] | Callable[[], float] | float, du0: list[float] | Callable[[], float] | float, ddu0: list[float] | Callable[[], float] | float) None#
class FEM.Core.CoreParabolic(geometry: Geometry, solver: Solver | None = None, verbose: bool = False, name='')#

Bases: CoreTransient

docstring for CoreParabolic

postProcess(t0, tf, steps) None#

Post process the solution

set_alpha(alpha)#
class FEM.Core.CoreTransient(geometry: Geometry, solver: Solver | None = None, verbose: bool = False, name='')#

Bases: Core

docstring for CoreTransient

apply_initial_condition(t0, dt)#
ensembling() None#

Ensembling of equation system. This method use the element gdl and the element matrices. The element matrices degrees of fredom must match the dimension of the element gdl. For m>1 variables per node, the gdl will be flattened. This ensure that the element matrices will always be a 2-D Numpy Array.

set_initial_condition(ic: list[float] | Callable[[], float] | float) None#

FEM.EDO1D module#

Solve a 1D 1 variable per node Finite Element problem

class FEM.EDO1D.EDO1D(geometry: Geometry, a: Callable, c: Callable, f: Callable, **kargs)#

Bases: Core

Create a 1D 1 variable per node Finite Element problem

The differential equation is:

\[a(x)\frac{d^2u}{dx^2}+c(x)u=f(x)\]
Parameters:
  • geometry (Geometry) – 1D Geometry of the problem. Use the Mesh.Lineal class

  • a (function) – Function a, if a is constant you can use a = lambda x: [value]

  • c (function) – Function c, if c is constant you can use c = lambda x: [value]

  • f (function) – Function f, if f is constant you can use f = lambda x: [value]

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model. Element matrices and forces are calculated with Gauss-Legendre quadrature. Point number depends of element discretization.

postProcess() None#

Generate graph of solution and solution derivative

FEM.Elasticity2D module#

2D Elasticity

class FEM.Elasticity2D.PlaneStrain(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStrain.<lambda>>, fy: ~typing.Callable = <function PlaneStrain.<lambda>>, **kargs)#

Bases: PlaneStress

Create a Plain Strain problem

Parameters:
  • geometry (Geometry) – 2D 2 variables per node geometry

  • E (int or float or list) – Young Moduli. If number, all element will have the same young moduli. If list, each position will be the element young moduli, so len(E) == len(self.elements)

  • v (int or float or list) – Poisson ratio. If number, all element will have the same Poisson ratio. If list, each position will be the element Poisson ratio, so len(v) == len(self.elements)

  • rho (int or float or list, optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (function, optional) – Function fx, if fx is constant you can use fx = lambda x: [value]. Defaults to lambda x:0.

  • fy (function, optional) – Function fy, if fy is constant you can use fy = lambda x: [value]. Defaults to lambda x:0.

class FEM.Elasticity2D.PlaneStrainSparse(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStrainSparse.<lambda>>, fy: ~typing.Callable = <function PlaneStrainSparse.<lambda>>, **kargs)#

Bases: PlaneStressSparse

Create a Plain Strain problem using sparse matrix

Parameters:
  • geometry (Geometry) – 2D 2 variables per node geometry with fast elements

  • E (int or float or list) – Young Moduli. If number, all element will have the same young moduli. If list, each position will be the element young moduli, so len(E) == len(self.elements)

  • v (int or float or list) – Poisson ratio. If number, all element will have the same Poisson ratio. If list, each position will be the element Poisson ratio, so len(v) == len(self.elements)

  • rho (int or float or list, optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (function, optional) – Function fx, if fx is constant you can use fx = lambda x: [value]. Defaults to lambda x:0.

  • fy (function, optional) – Function fy, if fy is constant you can use fy = lambda x: [value]. Defaults to lambda x:0.

class FEM.Elasticity2D.PlaneStress(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], t: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStress.<lambda>>, fy: ~typing.Callable = <function PlaneStress.<lambda>>, **kargs)#

Bases: PlaneStressOrthotropic

Create a Plain Stress problem

Parameters:
  • geometry (Geometry) – 2D 2 variables per node geometry

  • E (int or float or list) – Young Moduli. If number, all element will have the same young moduli. If list, each position will be the element young moduli, so len(E) == len(self.elements)

  • v (int or float or list) – Poisson ratio. If number, all element will have the same Poisson ratio. If list, each position will be the element Poisson ratio, so len(v) == len(self.elements)

  • t (int or float or list) – Element thickness. If number, all element will have the same thickness. If list, each position will be the element thickness, so len(t) == len(self.elements)

  • rho (int or float or list, optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (function, optional) – Function fx, if fx is constant you can use fx = lambda x: [value]. Defaults to lambda x:0.

  • fy (function, optional) – Function fy, if fy is constant you can use fy = lambda x: [value]. Defaults to lambda x:0.

class FEM.Elasticity2D.PlaneStressNonLocalSparse(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], t: ~typing.Tuple[float, list], l: float, z1: float, Lr: float, af: ~typing.Callable, rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStressNonLocalSparse.<lambda>>, fy: ~typing.Callable = <function PlaneStressNonLocalSparse.<lambda>>, notCalculateNonLocal=True, **kargs)#

Bases: PlaneStressSparse

Create a Plain Stress nonlocal problem using sparse matrices and the Pisano 2006 formulation.

Parameters:
  • geometry (Geometry) – 2D 2 variables per node geometry

  • E (int or float or list) – Young Moduli. If number, all element will have the same young moduli. If list, each position will be the element young moduli, so len(E) == len(self.elements)

  • v (int or float or list) – Poisson ratio. If number, all element will have the same Poisson ratio. If list, each position will be the element Poisson ratio, so len(v) == len(self.elements)

  • t (int or float or list) – Element thickness. If number, all element will have the same thickness. If list, each position will be the element thickness, so len(t) == len(self.elements)

  • l (float) – Internal lenght

  • z1 (float) – z1 factor

  • Lr (float) – Influence distance Lr

  • af (Callable) – Atenuation function

  • rho (int or float or list, optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (function, optional) – Function fx, if fx is constant you can use fx = lambda x: [value]. Defaults to lambda x:0.

  • fy (function, optional) – Function fy, if fy is constant you can use fy = lambda x: [value]. Defaults to lambda x:0.

elementMatrices() None#

Calculate the elements matrices

elementMatrix(ee: Element) None#

Calculates a single element local and nonlocal matrices

Parameters:

ee (Element) – Element to be calculated

ensembling() None#

Creation of the system sparse matrix. Force vector is ensembled in integration method

postProcess(mult: float = 1000, gs=None, levels=1000, **kargs) None#

Generate the stress surfaces and displacement fields for the geometry

Parameters:
  • mult (int, optional) – Factor for displacements. Defaults to 1000.

  • gs (list, optional) – List of 4 gridSpec matplotlib objects. Defaults to None.

profile(p0: list, p1: list, n: float = 100, plot=True) None#

Generate a profile between selected points

Parameters:
  • p0 (list) – start point of the profile [x0,y0]

  • p1 (list) – end point of the profile [xf,yf]

  • n (int, optional) – NUmber of samples for graph. Defaults to 100.

class FEM.Elasticity2D.PlaneStressNonLocalSparseNonHomogeneous(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], t: ~typing.Tuple[float, list], l: float, alpha: float, Lr: float, af: ~typing.Callable, rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStressNonLocalSparseNonHomogeneous.<lambda>>, fy: ~typing.Callable = <function PlaneStressNonLocalSparseNonHomogeneous.<lambda>>, **kargs)#

Bases: PlaneStressSparse

elementMatrices() None#

Calculate the elements matrices

elementMatrix(ee: Element) None#

Calculates a single element local and nonlocal matrices

Parameters:

ee (int) – Element to be calculated

ensembling() None#

Creation of the system sparse matrix. Force vector is ensembled in integration method

postProcess(mult: float = 1000, gs=None, levels=1000, **kargs) None#

Generate the stress surfaces and displacement fields for the geometry

Parameters:
  • mult (int, optional) – Factor for displacements. Defaults to 1000.

  • gs (list, optional) – List of 4 gridSpec matplotlib objects. Defaults to None.

profile(p0: list, p1: list, n: float = 100) None#

Generate a profile between selected points

Parameters:
  • p0 (list) – start point of the profile [x0,y0]

  • p1 (list) – end point of the profile [xf,yf]

  • n (int, optional) – NUmber of samples for graph. Defaults to 100.

class FEM.Elasticity2D.PlaneStressOrthotropic(geometry: ~FEM.Geometry.Geometry.Geometry, E1: ~typing.Tuple[float, list], E2: ~typing.Tuple[float, list], G12: ~typing.Tuple[float, list], v12: ~typing.Tuple[float, list], t: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStressOrthotropic.<lambda>>, fy: ~typing.Callable = <function PlaneStressOrthotropic.<lambda>>, **kargs)#

Bases: Core

Creates a plane stress problem with orthotropic formulation

Parameters:
  • geometry (Geometry) – Input geometry

  • E1 (Tuple[float, list]) – Young moduli in direction 1 (x)

  • E2 (Tuple[float, list]) – Young moduli in direction 2 (y)

  • G12 (Tuple[float, list]) – Shear moduli

  • v12 (Tuple[float, list]) – Poisson moduli

  • t (Tuple[float, list]) – Thickness

  • rho (Tuple[float, list], optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (Callable, optional) – Force in x direction. Defaults to lambdax:0.

  • fy (Callable, optional) – Force in y direction. Defaults to lambdax:0.

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

giveStressPoint(X: ndarray) Tuple[tuple, None]#

Calculates the stress in a given set of points.

Parameters:

X (np.ndarray) – Points to calculate the Stress. 2D Matrix. with 2 rows. First row is an array of 1 column with X coordinate. Second row is an array of 1 column with Y coordinate

Returns:

Tuple of stress (\(\sigma_x,\sigma_y,\sigma_{xy}\)) if X,Y exists in domain.

Return type:

tuple or None

postProcess(mult: float = 1000, gs=None, levels=1000, **kargs) None#

Generate the stress surfaces and displacement fields for the geometry

Parameters:
  • mult (int, optional) – Factor for displacements. Defaults to 1000.

  • gs (list, optional) – List of 4 gridSpec matplotlib objects. Defaults to None.

profile(p0: list, p1: list, n: float = 100) None#

Generate a profile between selected points

Parameters:
  • p0 (list) – start point of the profile [x0,y0]

  • p1 (list) – end point of the profile [xf,yf]

  • n (int, optional) – NUmber of samples for graph. Defaults to 100.

class FEM.Elasticity2D.PlaneStressOrthotropicSparse(geometry: ~FEM.Geometry.Geometry.Geometry, E1: ~typing.Tuple[float, list], E2: ~typing.Tuple[float, list], G12: ~typing.Tuple[float, list], v12: ~typing.Tuple[float, list], t: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStressOrthotropicSparse.<lambda>>, fy: ~typing.Callable = <function PlaneStressOrthotropicSparse.<lambda>>, **kargs)#

Bases: PlaneStressOrthotropic

Creates a plane stress problem with orthotropic formulation and sparce matrices

Parameters:
  • geometry (Geometry) – Input geometry

  • E1 (Tuple[float, list]) – Young moduli in direction 1 (x)

  • E2 (Tuple[float, list]) – Young moduli in direction 2 (y)

  • G12 (Tuple[float, list]) – Shear moduli

  • v12 (Tuple[float, list]) – Poisson moduli

  • t (Tuple[float, list]) – Thickness

  • rho (Tuple[float, list], optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (Callable, optional) – Force in x direction. Defaults to lambdax:0.

  • fy (Callable, optional) – Force in y direction. Defaults to lambdax:0.

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

ensembling() None#

Creation of the system sparse matrix. Force vector is ensembled in integration method

class FEM.Elasticity2D.PlaneStressSparse(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], t: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list] | None = None, fx: ~typing.Callable = <function PlaneStressSparse.<lambda>>, fy: ~typing.Callable = <function PlaneStressSparse.<lambda>>, **kargs)#

Bases: PlaneStressOrthotropicSparse

Create a Plain Stress problem using sparse matrices

Parameters:
  • geometry (Geometry) – 2D 2 variables per node geometry

  • E (int or float or list) – Young Moduli. If number, all element will have the same young moduli. If list, each position will be the element young moduli, so len(E) == len(self.elements)

  • v (int or float or list) – Poisson ratio. If number, all element will have the same Poisson ratio. If list, each position will be the element Poisson ratio, so len(v) == len(self.elements)

  • t (int or float or list) – Element thickness. If number, all element will have the same thickness. If list, each position will be the element thickness, so len(t) == len(self.elements)

  • rho (int or float or list, optional) – Density. If not given, mass matrix will not be calculated. Defaults to None.

  • fx (function, optional) – Function fx, if fx is constant you can use fx = lambda x: [value]. Defaults to lambda x:0.

  • fy (function, optional) – Function fy, if fy is constant you can use fy = lambda x: [value]. Defaults to lambda x:0.

FEM.Elasticity3D module#

3D Isotropic elasticity

class FEM.Elasticity3D.Elasticity(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list], fx: ~typing.Callable = <function Elasticity.<lambda>>, fy: ~typing.Callable = <function Elasticity.<lambda>>, fz: ~typing.Callable = <function Elasticity.<lambda>>, **kargs)#

Bases: Core

Creates a 3D Elasticity problem

Parameters:
  • geometry (Geometry) – Input 3D geometry

  • E (Tuple[float, list]) – Young Moduli

  • v (Tuple[float, list]) – Poisson coeficient

  • rho (Tuple[float, list]) – Density.

  • fx (Callable, optional) – Force in x direction. Defaults to lambdax:0.

  • fy (Callable, optional) – Force in y direction. Defaults to lambdax:0.

  • fz (Callable, optional) – Force in z direction. Defaults to lambdax:0.

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

ensembling() None#

Creation of the system sparse matrix. Force vector is ensembled in integration method

postProcess(**kargs) None#

Calculate stress and strain for each element gauss point and vertice. The results are stored in each element sigmas and epsilons properties.

profile(region: list[float], n: float = 10) None#

Creates a profile in a given region coordinates

Parameters:
  • region (list[float]) – List of region coordinates (square region 2D)

  • n (float, optional) – Number of points. Defaults to 10.

Returns:

Coordinates, displacements and second variable solution

Return type:

np.ndarray

class FEM.Elasticity3D.ElasticityFromTensor(geometry: ~FEM.Geometry.Geometry.Geometry, C: ~numpy.ndarray, rho: ~typing.Tuple[float, list], fx: ~typing.Callable = <function ElasticityFromTensor.<lambda>>, fy: ~typing.Callable = <function ElasticityFromTensor.<lambda>>, fz: ~typing.Callable = <function ElasticityFromTensor.<lambda>>, **kargs)#

Bases: Elasticity

class FEM.Elasticity3D.NonLocalElasticity(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list], l: float, z1: float, Lr: float, af: ~typing.Callable, fx: ~typing.Callable = <function NonLocalElasticity.<lambda>>, fy: ~typing.Callable = <function NonLocalElasticity.<lambda>>, fz: ~typing.Callable = <function NonLocalElasticity.<lambda>>, **kargs)#

Bases: Elasticity

Creates a 3D Elasticity problem

Parameters:
  • geometry (Geometry) – Input 3D geometry

  • E (Tuple[float, list]) – Young Moduli

  • v (Tuple[float, list]) – Poisson coeficient

  • rho (Tuple[float, list]) – Density.

  • l (float) – Internal lenght

  • z1 (float) – Z1 factor

  • Lr (float) – Influence distance

  • af (Callable) – Atenuation function

  • fx (Callable, optional) – Force in x direction. Defaults to lambdax:0.

  • fy (Callable, optional) – Force in y direction. Defaults to lambdax:0.

  • fz (Callable, optional) – Force in z direction. Defaults to lambdax:0.

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

ensembling() None#

Creation of the system sparse matrix. Force vector is ensembled in integration method

class FEM.Elasticity3D.NonLocalElasticityFromTensor(geometry: ~FEM.Geometry.Geometry.Geometry, C: ~numpy.ndarray, rho: ~typing.Tuple[float, list], l: float, z1: float, Lr: float, af: ~typing.Callable, fx: ~typing.Callable = <function NonLocalElasticityFromTensor.<lambda>>, fy: ~typing.Callable = <function NonLocalElasticityFromTensor.<lambda>>, fz: ~typing.Callable = <function NonLocalElasticityFromTensor.<lambda>>, **kargs)#

Bases: NonLocalElasticity

class FEM.Elasticity3D.NonLocalElasticityLegacy(geometry: ~FEM.Geometry.Geometry.Geometry, E: ~typing.Tuple[float, list], v: ~typing.Tuple[float, list], rho: ~typing.Tuple[float, list], l: float, z1: float, Lr: float, af: ~typing.Callable, fx: ~typing.Callable = <function NonLocalElasticityLegacy.<lambda>>, fy: ~typing.Callable = <function NonLocalElasticityLegacy.<lambda>>, fz: ~typing.Callable = <function NonLocalElasticityLegacy.<lambda>>, **kargs)#

Bases: Elasticity

Creates a 3D Elasticity problem

Parameters:
  • geometry (Geometry) – Input 3D geometry

  • E (Tuple[float, list]) – Young Moduli

  • v (Tuple[float, list]) – Poisson coeficient

  • rho (Tuple[float, list]) – Density.

  • l (float) – Internal lenght

  • z1 (float) – Z1 factor

  • Lr (float) – Influence distance

  • af (Callable) – Atenuation function

  • fx (Callable, optional) – Force in x direction. Defaults to lambdax:0.

  • fy (Callable, optional) – Force in y direction. Defaults to lambdax:0.

  • fz (Callable, optional) – Force in z direction. Defaults to lambdax:0.

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

ensembling() None#

Creation of the system sparse matrix. Force vector is ensembled in integration method

FEM.EulerBernoulliBeam module#

Euler Bernoulli Beam implementation

class FEM.EulerBernoulliBeam.EulerBernoulliBeam(geometry: Geometry, EI: float, cf: float = 0.0, f: float = 0.0, **kargs)#

Bases: Core

Creates a Euler Bernoulli beam problem

Parameters:
  • geometry (Geometry) – 1D 2 variables per node problem geometry. Geometry must have Euler Bernoulli elements.

  • EI (float) – Young’s moduli multiplied by second moment of area (inertia).

  • cf (float, optional) – Soil coeficient. Defaults to 0.

  • f (float or int or function, optional) – Force function applied to the beam. Defaults to 0.0

elementMatrices() None#

Calculate the element matrices usign Guass Legendre quadrature.

postProcess(plot=True) None#

Post process the solution. Shows graphs of displacement, rotation, shear and moment.

class FEM.EulerBernoulliBeam.EulerBernoulliBeamNonLineal(geometry: Geometry, EI: float, EA: float, fx: float = 0.0, fy: float = 0.0)#

Bases: Core

Creates a Euler Bernoulli beam problem

Parameters:
  • geometry (Geometry) – 1D 2 variables per node problem geometry. Geometry must have Euler Bernoulli elements.

  • EI (float) – Young’s moduli multiplied by second moment of area (inertia).

  • EA (float) – Young’s moduli multiplied by area.

  • cf (float, optional) – Soil coeficient. Defaults to 0.

  • fx (float or int or function, optional) – Force function applied in the x direction to the beam. Defaults to 0.0

  • fy (float or int or function, optional) – Force function applied in the y direction to the beam. Defaults to 0.0

elementMatrices() None#

Calculate the element matrices usign Guass Legendre quadrature.

postProcess(plot=True) None#

Post process the solution. Shows graphs of displacement, rotation, shear and moment.

FEM.FEMLogger module#

class FEM.FEMLogger.FEMLogger(ad: str = '')#

Bases: object

Creation of a Logger for FEM purposes. Based on Python Logger by Fonic <https://github.com/fonic>

end_timer()#

Ends the sesion time

setup_logging(console_log_output='stdout', console_log_level='warning', console_log_color=True, logfile_file=None, logfile_log_level='debug', logfile_log_color=False, log_line_template='%(color_on)s[%(levelname)-8s] %(message)s%(color_off)s')#

Set the logger

Parameters:
  • console_log_output (str, optional) – . Defaults to “stdout”.

  • console_log_level (str, optional) – . Defaults to “warning”.

  • console_log_color (bool, optional) – . Defaults to True.

  • logfile_file (optional) – . Defaults to None.

  • logfile_log_level (str, optional) – . Defaults to “debug”.

  • logfile_log_color (bool, optional) – . Defaults to False.

  • log_line_template (str, optional) – . Defaults to “%(color_on)s[%(levelname)-8s] %(message)s%(color_off)s”.

class FEM.FEMLogger.LogFormatter(color, *args, **kwargs)#

Bases: Formatter

Creates a Logging Formatter

Parameters:

color (color) – Color

COLOR_CODES = {10: '\x1b[1;30m', 20: '\x1b[0;37m', 30: '\x1b[1;33m', 40: '\x1b[1;31m', 50: '\x1b[1;35m'}#
RESET_CODE = '\x1b[0m'#
format(record, *args, **kwargs)#

Formats arecord

Parameters:

record (record) – Record to be formatted

class FEM.FEMLogger.TimeFilter(name='')#

Bases: Filter

filter(record)#

Determine if the specified record is to be logged.

Returns True if the record should be logged, or False otherwise. If deemed appropriate, the record may be modified in-place.

FEM.Heat1D module#

1D Stady state heat with convective border conditions

class FEM.Heat1D.Heat1D(geometry: Geometry, A: float, P: float, ku: float, beta: float, Ta: float, q: float = 0.0, **kargs)#

Bases: Core

Creates a 1D Stady state heat problem. Convective border conditions can be applied

The differential equation is:

\[-\frac{d}{dx}\left(Ak\frac{dT}{dx}\right)+\beta P(T-T_{\infty})=q\]
Parameters:
  • geometry (Geometry) – Input 1D Geometry. 1 variable per node

  • A (float or list) – Section area. If float, all elements will have the same area. If list, the i-element will have the i-area

  • P (float or list) – Section perimeter. If float, all elements will have the same perimeter. If list, the i-element will have the i-perimeter

  • ku (float or list) – Conductivity. If float, all elements will have the same conductivity. If list, the i-element will have the i-conductivity

  • beta (float or list) – Transfer coeficient. If float, all elements will have the same transfer coeficient. If list, the i-element will have the i-transfer coeficient

  • Ta (float) – Ambient temperature

  • q (float or list, optional) – Internal heat generation rate. If float, all elements will have the same internal heat generation rate coeficient. If list, the i-element will have the i-internal heat generation rate. Defaults to 0.0.

borderConditions() None#

Assign border conditions to the system. The border conditios are assigned in this order:

  1. Natural border conditions

  2. Essential border conditions

This ensures that in a node with 2 border conditions the essential border conditions will be applied.

defineConvectiveBoderConditions(node: int, value: float = 0) None#

Add a convective border condition. The value is: \(kA\frac{dT}{dx}+\beta A(T-T_{\infty})=value\)

Parameters:
  • node (int) – Node where the above border condition is applied

  • value (float, optional) – Defined below. Defaults to 0.

elementMatrices() None#

Calculate the element matrices using Gauss Legendre quadrature.

postProcess() None#

Generate graph of solution and solution derivative

set_convective_conditions()#
class FEM.Heat1D.Heat1DTransient(geometry: Geometry, A: float, P: float, ku: float, beta: float, Ta: float, q: float = 0.0, **kargs)#

Bases: Heat1D, CoreParabolic

Creates a 1D Stady state heat problem. Convective border conditions can be applied

The differential equation is:

\[\frac{\partial T}{\partial t}-\frac{\partial}{\partial x}\left(Ak\frac{\partial T}{\partial x}\right)+\beta P(T-T_{\infty})=q\]
Parameters:
  • geometry (Geometry) – Input 1D Geometry. 1 variable per node

  • A (float or list) – Section area. If float, all elements will have the same area. If list, the i-element will have the i-area

  • P (float or list) – Section perimeter. If float, all elements will have the same perimeter. If list, the i-element will have the i-perimeter

  • k (float or list) – Conductivity. If float, all elements will have the same conductivity. If list, the i-element will have the i-conductivity

  • beta (float or list) – Transfer coeficient. If float, all elements will have the same transfer coeficient. If list, the i-element will have the i-transfer coeficient

  • Ta (float) – Ambient temperature (also called T∞)

  • q (float or list, optional) – Internal heat generation rate. If float, all elements will have the same internal heat generation rate coeficient. If list, the i-element will have the i-internal heat generation rate. Defaults to 0.0.

elementMatrices() None#

Calculate the element matrices using Gauss Legendre quadrature.

postProcess(node=None, t0=None, tf=None, steps=None, dt=None, ax=None) None#

Post process the solution and steps

FEM.Heat2D module#

2D Stady state heat with convective border conditions

class FEM.Heat2D.Heat2D(geometry: Geometry, kx: Tuple[float, list], ky: Tuple[float, list], f: Callable | None = None, **kargs)#

Bases: Core

Creates a Heat2D problem with convective borders

The differential equation is:

\[-\frac{\partial}{\partial x}\left(k_x\frac{\partial T}{\partial x}\right) - \frac{\partial}{\partial y}\left(k_y\frac{\partial T}{\partial y}\right)=f(x,y)\]

With convective border conditions:

\[k_x\frac{\partial T}{\partial x}n_x+k_y\frac{\partial T}{\partial y}n_y+\beta (T-T_\infty)=\hat{q_n}\]
Parameters:
  • geometry (Geometry) – Input 1 variable per node geometry

  • kx (Tuple[float, list]) – Heat transfer coeficient in x direction. If number, all element will have the same coefficient. If list, each position will be the element coefficient, so len(kx) == len(self.elements)

  • ky (Tuple[float, list]) – Heat transfer coeficient in y direction. If number, all element will have the same coefficient. If list, each position will be the element coefficient, so len(kx) == len(self.elements)

  • f (Callable, optional) – Internal heat generation function. Defaults to None.

defineConvectiveBoderConditions(region: int, beta: float = 0, Ta: float = 0) None#

Define convective borders

Parameters:
  • region (int) – region in wich load will be applied

  • beta (float, optional) – Convective coeficient \(\beta\) . Defaults to 0.

  • Ta (float, optional) – Ambient temperature in convective border. Defaults to 0.

elementMatrices() None#

Calculate the element matrices using Gauss Legendre quadrature.

postProcess(levels=1000) None#

Generate the temperature surface for the geometry

FEM.NonLinealExample module#

Solve a 1D 1 variable per node non lineal Finite Element problem

class FEM.NonLinealExample.NonLinealSimpleEquation(geometry: Geometry, a: Callable, f: Callable, **kargs)#

Bases: Core

Creates a nonlineal 1D equation with the form:

\[-\frac{d}{dx}\left(a(x)u\frac{du}{dx}\right)=f(x)\]
Parameters:
  • geometry (Geometry) – Input lineal geometry

  • a (Callable) – Function a

  • f (Callable) – Function f

elementMatrices() None#

Calculate the element matrices usign Reddy’s non lineal finite element model. Element matrices and forces are calculated with Gauss-Legendre quadrature. Point number depends of element discretization.

postProcess() None#

Generate graph of solution and solution derivative

FEM.Poisson2D module#

Solve the Poisson equation for a 2D domain.

class FEM.Poisson2D.Poisson2D(geometry: Geometry, phi: float)#

Bases: Core

Create a Poisson2D finite element problem

The differential equation is:

\[-\nabla^2\Psi=\theta\]
Parameters:
  • geometry (Geometry) – 2D 1 variable per node geometry

  • phi (float) – Value

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

postProcess(levels=30, derivatives=True) None#

Create graphs for solution and derivatives.

FEM.Torsion2D module#

Solve the torsion Poisson equation for a 2D domain.

class FEM.Torsion2D.Torsion2D(geometry: Geometry, G: float, phi: float, **kargs)#

Bases: Core

Create a torsional finite element problem

The differential equation is:

\[-\frac{-1}{G}\nabla^2\Psi=2\theta\]

Where:

\[\begin{split}\sigma_{xz}\frac{\partial\Psi}{\partial y} \\ \sigma_{yz}=-\frac{\partial\Psi}{\partial x}\end{split}\]

With \(\Psi=0\) on the boundary.

Parameters:
  • geometry (Geometry) – 2D 1 variable per node geometry

  • G (float) – Shear moduli of elements

  • phi (float) – Rotation angle in radians

elementMatrices() None#

Calculate the element matrices usign Reddy’s (2005) finite element model

postProcess(levels=30, derivatives=True) None#

Create graphs for stress function and derivatives.