Finite Element Solvers

1D Poisson Equation

class FEMpy.Solvers.Poisson1D(mesh, fe_trial_basis, fe_test_basis, boundary_conditions)[source]

Solves a one-dimensional Poisson equation.

Uses finite element methods to solve a Poisson differential equation of the form

\[- \frac{{\rm d}}{{\rm d} x}\left(c(x) \frac{{\rm d}}{{\rm d} x} u(x) \right) = f(x); a \leq x \leq b.\]

with a combination of Dirichlet, Neumann, or Robin boundary conditions.

Warning

Both end point boundary conditions cannot be Neumann as this may result in a loss of uniqueness of the solution.

Parameters:

Examples

>>> import numpy as np
>>> from FEMpy import Interval1D, IntervalBasis1D, BoundaryConditions, Poisson1D
>>> mesh = Interval1D(0, 1, h=1/2, basis_type='linear')
>>> basis = IntervalBasis1D('linear')
>>> dirichlet_funct = lambda x: 0 if x == 0 else np.cos(1)
>>> bcs = BoundaryConditions(mesh, ('dirichlet', 'dirichlet'), dirichlet_fun=dirichlet_funct)
>>> coefficient_funct = lambda x: np.exp(x)
>>> source_funct = lambda x: -np.exp(x) * (np.cos(x) - 2*np.sin(x) - x*np.cos(x) - x*np.sin(x))
>>> poisson_eq = Poisson1D(mesh, basis, basis, bcs)
>>> poisson_eq.solve(coefficient_funct, source_funct)
array([0.        , 0.44814801, 0.54030231])
fe_solution(x, local_sol, vertices, derivative_order)[source]

Defines the functional solution piecewise on the finite on the finite elements.

Uses the solution vector and the basis function to define a piecewise continuous solution over the element.

Parameters:
  • x (float or array_like) – A value or array of points to evaluate the function on.
  • local_sol (array_like) – Finite element solution node vector local to the element En.
  • vertices (array_like) – Global node coordinates for the mesh element En.
  • derivative_order (int) – The derivative order to take the basis function to.
Returns:

Solution at all points in x in element.

Return type:

float

h1_seminorm_error(diff_exact_sol)[source]

The H1 semi-norm error of the finite element solution compared against the given analyatical solution.

Parameters:diff_exact_sol (function) – The first derivative of the analytical solution to the Poisson equation.
Returns:The H1 semi-norm error of the finite element solution over the domain evaluated element-wise.
Return type:float
l2_error(exact_sol)[source]

The L2 norm error of the finite element solution compared against the given analytical solution.

Parameters:exact_sol (function) – The analytical solution to the Poisson equation.
Returns:The L2 norm error of the finite element solution over the domain evaluated element-wise.
Return type:float
l_inf_error(exact_sol)[source]

Computes the L-infinity norm error.

Computes the L-infinity norm error using the exact solution and the finite element function fe_solution.

Parameters:exact_sol (function) – The analytical solution to compare the finite element solution against.
Returns:The L-infinity norm error of the finite element solution over the domain evaluated element-wise.
Return type:float
solve(coeff_fun, source_fun)[source]

Method that performs the finte element solution algorithm.

Calls the assembly functions FEMpy.Assemblers.assemble_matrix and FEMpy.Assemblers.assemble_vector to create the stiffness matrix and load vector respectively. Then, applies the boundary condition treatments to the matrix and vector. Finally, solves the linear system \(A\mathbf{x} = \mathbf{b}.\)

Parameters:
  • coeff_fun (function) – Function name of the coefficient function c`(`x) in the Poisson equation.
  • source_fun (function) – The nonhomogeneous source function f`(`x) of the Poisson equation.
Returns:

The nodal solution vector.

Return type:

ndarray

2D Poisson Equation

class FEMpy.Solvers.Poisson2D(mesh, fe_trial_basis, fe_test_basis, boundary_conditions)[source]

Solves a two-dimensional Poisson equation.

Uses finite element methods to solve a Poisson differential equation of the form

\[-\nabla\left(c(\mathbf{x}) \cdot \nabla u(\mathbf{x}) \right) = f(\mathbf{x}); \mathbf{x} \in \Omega\]

with a combination of Dirichlet, Neumann, or Robin boundary conditions.

Warning

All edge boundary conditions cannot be Neumann as this may result in a loss of uniqueness of the solution.

Parameters:

Examples

>>> import numpy as np
>>> from FEMpy import TriangularMesh2D, TriangularBasis2D, BoundaryConditions2D, Poisson2D
>>> left, right, bottom, top = -1, 1, -1, 1
>>> h = 1
>>> def dirichlet_funct(coord):
>>>     x, y = coord
>>>     if x == -1:
>>>         return np.exp(-1 + y)
>>>     elif x == 1:
>>>         return np.exp(1 + y)
>>>     elif y == 1:
>>>         return np.exp(x + 1)
>>>     elif y == -1:
>>>         return np.exp(x - 1)
>>> coeff_funct = lambda coord: 1
>>> source_funct = lambda coord: -2 * np.exp(coord[0] + coord[1])
>>> mesh = TriangularMesh2D(left, right, bottom, top, h, h, 'linear')
>>> basis = TriangularBasis2D('linear')
>>> boundary_node_types = ['dirichlet'] * mesh.boundary_nodes.shape[1]
>>> boundary_edge_types = ['dirichlet'] * (mesh.boundary_edges.shape[1]-1)
>>> bcs = BoundaryConditions2D(mesh, boundary_node_types, boundary_edge_types, dirichlet_fun=dirichlet_funct)
>>> poisson_eq = Poisson2D(mesh, basis, basis, bcs)
>>> poisson_eq.solve(coeff_funct, source_funct)
array([0.13533528, 0.36787944, 1., 0.36787944, 1.,  2.71828183, 1., 2.71828183, 7.3890561])
fe_solution(coords, local_sol, vertices, derivative_order)[source]

Defines the functional solution piecewise on the finite on the finite elements.

Uses the solution vector and the basis function to define a piecewise continuous solution over the element.

Parameters:
  • coords (float or array_like) – A value or array of points to evaluate the function on.
  • local_sol (array_like) – Finite element solution node vector local to the element En.
  • vertices (array_like) – Global node coordinates for the mesh element En.
  • derivative_order (tuple of int) – The derivative orders in the x- and y-directions to take the basis function to.
Returns:

Solution at all points in coords in element.

Return type:

float

h1_seminorm_error(diff_exact_sol)[source]

The H1 semi-norm error of the finite element solution compared against the given analyatical solution.

Parameters:diff_exact_sol (tuple of function) – A tuple of first derivatives in the x- and the y- directions the analytical solution to the Poisson equation respectively.
Returns:The full H1 semi-norm error of the finite element solution over the domain evaluated element-wise.
Return type:float
l2_error(exact_sol)[source]

The L2 norm error of the finite element solution compared against the given analytical solution.

Parameters:exact_sol (function) – The analytical solution to the Poisson equation.
Returns:The L2 norm error of the finite element solution over the domain evaluated element-wise.
Return type:float
l_inf_error(exact_sol)[source]

Computes the L-infinity norm error.

Computes the L-infinity norm error using the exact solution and the finite element function fe_solution.

Parameters:exact_sol (function) – The analytical solution to compare the finite element solution against.
Returns:The L-infinity norm error of the finite element solution over the domain evaluated element-wise.
Return type:float
solve(coeff_fun, source_fun)[source]

Method that performs the finte element solution algorithm.

Calls the assembly functions FEMpy.Assemblers.assemble_matrix and FEMpy.Assemblers.assemble_vector to create the stiffness matrix and load vector respectively. Then, applies the boundary condition treatments to the matrix and vector. Finally, solves the linear system \(A\mathbf{x} = \mathbf{b}.\)

Parameters:
  • coeff_fun (function) – Function name of the coefficient function c`(`x, y) in the Poisson equation.
  • source_fun (function) – The nonhomogeneous source function f`(`x, y) of the Poisson equation.