implementations.problem_classes.HeatEquation_1D_FEniCS_matrix_forced module

class fenics_heat(c_nvars=128, t0=0.0, family='CG', order=4, refinements=1, nu=0.1, c=0.0)[source]

Bases: ptype

Example implementing the forced one-dimensional heat equation with Dirichlet boundary conditions

\[\frac{d u}{d t} = \nu \frac{d^2 u}{d x^2} + f\]

for \(x \in \Omega:=[0,1]\), where the forcing term \(f\) is defined by

\[f(x, t) = -\sin(\pi x) (\sin(t) - \nu \pi^2 \cos(t)).\]

For initial conditions with constant c and

\[u(x, 0) = \sin(\pi x) + c\]

the exact solution of the problem is given by

\[u(x, t) = \sin(\pi x)\cos(t) + c.\]

In this class the problem is implemented in the way that the spatial part is solved using FEniCS [1]_. Hence, the problem is reformulated to the weak formulation

The part containing the forcing term is treated explicitly, where it is interpolated in the function space. The other part will be treated in an implicit way.

Parameters:
  • c_nvars (int, optional) – Spatial resolution, i.e., numbers of degrees of freedom in space.

  • t0 (float, optional) – Starting time.

  • family (str, optional) – Indicates the family of elements used to create the function space for the trail and test functions. The default is 'CG', which are the class of Continuous Galerkin, a synonym for the Lagrange family of elements, see [2]_.

  • order (int, optional) – Defines the order of the elements in the function space.

  • refinements (int, optional) – Denotes the refinement of the mesh. refinements=2 refines the mesh by factor \(2\).

  • nu (float, optional) – Diffusion coefficient \(\nu\).

  • c (float, optional) – Constant for the Dirichlet boundary condition :math: c

V

Defines the function space of the trial and test functions.

Type:

FunctionSpace

M

Denotes the expression \(\int_\Omega u_t v\,dx\).

Type:

scalar, vector, matrix or higher rank tensor

K

Denotes the expression \(- \nu \int_\Omega \nabla u \nabla v\,dx\).

Type:

scalar, vector, matrix or higher rank tensor

g

The forcing term \(f\) in the heat equation.

Type:

Expression

bc

Denotes the Dirichlet boundary conditions.

Type:

DirichletBC

References

apply_mass_matrix(u)[source]

Routine to apply mass matrix.

Parameters:

u (dtype_u) – Current values of the numerical solution.

Returns:

me – The product \(M \vec{u}\).

Return type:

dtype_u

dtype_f

alias of rhs_fenics_mesh

dtype_u

alias of fenics_mesh

eval_f(u, t)[source]

Routine to evaluate both parts of the right-hand side of the problem.

Parameters:
  • u (dtype_u) – Current values of the numerical solution.

  • t (float) – Current time at which the numerical solution is computed.

Returns:

f – The right-hand side divided into two parts.

Return type:

dtype_f

solve_system(rhs, factor, u0, t)[source]

Dolfin’s linear solver for \((M - factor \cdot A) \vec{u} = \vec{rhs}\).

Parameters:
  • rhs (dtype_f) – Right-hand side for the nonlinear system.

  • factor (float) – Abbrev. for the node-to-node stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver (not used here so far).

  • t (float) – Current time.

Returns:

u – Solution.

Return type:

dtype_u

u_exact(t)[source]

Routine to compute the exact solution at time \(t\).

Parameters:

t (float) – Time of the exact solution.

Returns:

me – Exact solution.

Return type:

dtype_u

class fenics_heat_mass(c_nvars=128, t0=0.0, family='CG', order=4, refinements=1, nu=0.1, c=0.0)[source]

Bases: fenics_heat

Example implementing the forced one-dimensional heat equation with Dirichlet boundary conditions

\[\frac{d u}{d t} = \nu \frac{d^2 u}{d x^2} + f\]

for \(x \in \Omega:=[0,1]\), where the forcing term \(f\) is defined by

\[f(x, t) = -\sin(\pi x) (\sin(t) - \nu \pi^2 \cos(t)).\]

For initial conditions with constant c and

\[u(x, 0) = \sin(\pi x) + c\]

the exact solution of the problem is given by

\[u(x, t) = \sin(\pi x)\cos(t) + c.\]

In this class the problem is implemented in the way that the spatial part is solved using FEniCS [1]_. Hence, the problem is reformulated to the weak formulation

The forcing term is treated explicitly, and is expressed via the mass matrix resulting from the left-hand side term \(\int_\Omega u_t v\,dx\), and the other part will be treated in an implicit way.

Parameters:
  • c_nvars (int, optional) – Spatial resolution, i.e., numbers of degrees of freedom in space.

  • t0 (float, optional) – Starting time.

  • family (str, optional) – Indicates the family of elements used to create the function space for the trail and test functions. The default is 'CG', which are the class of Continuous Galerkin, a synonym for the Lagrange family of elements, see [2]_.

  • order (int, optional) – Defines the order of the elements in the function space.

  • refinements (int, optional) – Denotes the refinement of the mesh. refinements=2 refines the mesh by factor \(2\).

  • nu (float, optional) – Diffusion coefficient \(\nu\).

V

Defines the function space of the trial and test functions.

Type:

FunctionSpace

M

Denotes the expression \(\int_\Omega u_t v\,dx\).

Type:

scalar, vector, matrix or higher rank tensor

K

Denotes the expression \(- \nu \int_\Omega \nabla u \nabla v\,dx\).

Type:

scalar, vector, matrix or higher rank tensor

g

The forcing term \(f\) in the heat equation.

Type:

Expression

bc

Denotes the Dirichlet boundary conditions.

Type:

DirichletBC

bc_hom

Denotes the homogeneous Dirichlet boundary conditions, potentially required for fixing the residual

Type:

DirichletBC

fix_bc_for_residual

flag to indicate that the residual requires special treatment due to boundary conditions

Type:

boolean

References

eval_f(u, t)[source]

Routine to evaluate both parts of the right-hand side.

Parameters:
  • u (dtype_u) – Current values of the numerical solution.

  • t (float) – Current time at which the numerical solution is computed.

Returns:

f – The right-hand side divided into two parts.

Return type:

dtype_f

fix_residual(res)[source]

Applies homogeneous Dirichlet boundary conditions to the residual

Parameters:

res (dtype_u) – Residual

solve_system(rhs, factor, u0, t)[source]

Dolfin’s linear solver for \((M - factor A) \vec{u} = \vec{rhs}\).

Parameters:
  • rhs (dtype_f) – Right-hand side for the nonlinear system.

  • factor (float) – Abbrev. for the node-to-node stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver (not used here so far).

  • t (float) – Current time.

Returns:

u – Solution.

Return type:

dtype_u

class fenics_heat_mass_timebc(c_nvars=128, t0=0.0, family='CG', order=4, refinements=1, nu=0.1, c=0.0)[source]

Bases: fenics_heat_mass

Example implementing the forced one-dimensional heat equation with time-dependent Dirichlet boundary conditions

\[\frac{d u}{d t} = \nu \frac{d^2 u}{d x^2} + f\]

for \(x \in \Omega:=[0,1]\), where the forcing term \(f\) is defined by

\[f(x, t) = -\cos(\pi x) (\sin(t) - \nu \pi^2 \cos(t)).\]

and the boundary conditions are given by

\[u(x, t) = \cos(\pi x)\cos(t).\]

The exact solution of the problem is given by

\[u(x, t) = \cos(\pi x)\cos(t) + c.\]

In this class the problem is implemented in the way that the spatial part is solved using FEniCS [1]_. Hence, the problem is reformulated to the weak formulation

The forcing term is treated explicitly, and is expressed via the mass matrix resulting from the left-hand side term \(\int_\Omega u_t v\,dx\), and the other part will be treated in an implicit way.

Parameters:
  • c_nvars (int, optional) – Spatial resolution, i.e., numbers of degrees of freedom in space.

  • t0 (float, optional) – Starting time.

  • family (str, optional) – Indicates the family of elements used to create the function space for the trail and test functions. The default is 'CG', which are the class of Continuous Galerkin, a synonym for the Lagrange family of elements, see [2]_.

  • order (int, optional) – Defines the order of the elements in the function space.

  • refinements (int, optional) – Denotes the refinement of the mesh. refinements=2 refines the mesh by factor \(2\).

  • nu (float, optional) – Diffusion coefficient \(\nu\).

V

Defines the function space of the trial and test functions.

Type:

FunctionSpace

M

Denotes the expression \(\int_\Omega u_t v\,dx\).

Type:

scalar, vector, matrix or higher rank tensor

K

Denotes the expression \(- \nu \int_\Omega \nabla u \nabla v\,dx\).

Type:

scalar, vector, matrix or higher rank tensor

g

The forcing term \(f\) in the heat equation.

Type:

Expression

bc

Denotes the time-dependent Dirichlet boundary conditions.

Type:

DirichletBC

bc_hom

Denotes the homogeneous Dirichlet boundary conditions, potentially required for fixing the residual

Type:

DirichletBC

fix_bc_for_residual

flag to indicate that the residual requires special treatment due to boundary conditions

Type:

boolean

References

solve_system(rhs, factor, u0, t)[source]

Dolfin’s linear solver for \((M - factor A) \vec{u} = \vec{rhs}\).

Parameters:
  • rhs (dtype_f) – Right-hand side for the nonlinear system.

  • factor (float) – Abbrev. for the node-to-node stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver (not used here so far).

  • t (float) – Current time.

Returns:

u – Solution.

Return type:

dtype_u

u_exact(t)[source]

Routine to compute the exact solution at time \(t\).

Parameters:

t (float) – Time of the exact solution.

Returns:

me – Exact solution.

Return type:

dtype_u