implementations.problem_classes.GrayScott_1D_FEniCS_implicit module

class fenics_grayscott(c_nvars=256, t0=0.0, family='CG', order=4, refinements=None, Du=1.0, Dv=0.01, A=0.09, B=0.086)[source]

Bases: ptype

The Gray-Scott system [1] describes a reaction-diffusion process of two substances \(u\) and \(v\), where they diffuse over time. During the reaction \(u\) is used up with overall decay rate \(B\), whereas \(v\) is produced with feed rate \(A\). \(D_u,\, D_v\) are the diffusion rates for \(u,\, v\). This process is described by the one-dimensional model using Dirichlet boundary conditions

\[\frac{\partial u}{\partial t} = D_u \Delta u - u v^2 + A (1 - u),\]
\[\frac{\partial v}{\partial t} = D_v \Delta v + u v^2 - B u\]

for \(x \in \Omega:=[0, 100]\). The weak formulation of the problem can be obtained by multiplying the system with a test function \(q\):

\[\int_\Omega u_t q\,dx = \int_\Omega D_u \Delta u q - u v^2 q + A (1 - u) q\,dx,\]
\[\int_\Omega v_t q\,dx = \int_\Omega D_v \Delta v q + u v^2 q - B u q\,dx,\]

The spatial solve of the weak formulation is realized by FEniCS [2].

Parameters:
  • c_nvars (int, optional) – Spatial resolution, i.e., number 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 [3].

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

  • refinements (list or tuple, optional) – Defines the refinement for the spatial grid. Needs to be a list or tuple, e.g. refinements=[2, 2] or refinements=(2, 2).

  • Du (float, optional) – Diffusion rate for \(u\).

  • Dv (float, optional) – Diffusion rate for \(v\).

  • A (float, optional) – Feed rate for \(v\).

  • B (float, optional) – Overall decay rate for \(u\).

V

Defines the function space of the trial and test functions.

Type:

FunctionSpace

w

Function for the right-hand side.

Type:

Function

w1

Split of w, part 1.

Type:

Function

w2

Split of w, part 2.

Type:

Function

F1

Weak form of right-hand side, first part.

Type:

scalar, vector, matrix or higher rank tensor

F2

Weak form of right-hand side, second part.

Type:

scalar, vector, matrix or higher rank tensor

F

Weak form of full right-hand side.

Type:

scalar, vector, matrix or higher rank tensor

M

Full mass matrix for both parts.

Type:

matrix

References

dtype_f

alias of fenics_mesh

dtype_u

alias of fenics_mesh

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

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:

me – Solution as mesh.

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 (only at \(t_0 = 0.0\)).

Return type:

dtype_u