implementations.problem_classes.GrayScott_MPIFFT module

class grayscott_imex_diffusion(Du=1.0, Dv=0.01, A=0.09, B=0.086, **kwargs)[source]

Bases: IMEX_Laplacian_MPIFFT

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\). Here, the process is described by the \(N\)-dimensional model

\[\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\]

in \(x \in \Omega:=[-L/2, L/2]^N\) with \(N=2,3\). Spatial discretization is done by using Fast Fourier transformation for solving the linear parts provided by mpi4py-fft [2]_, see also https://mpi4py-fft.readthedocs.io/en/latest/.

This class implements the problem for semi-explicit time-stepping (diffusion is treated implicitly, and reaction is computed in explicit fashion).

Parameters:
  • nvars (tuple of int, optional) – Spatial resolution, i.e., number of degrees of freedom in space. Should be a tuple, e.g. nvars=(127, 127).

  • 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\).

  • spectral (bool, optional) – If True, the solution is computed in spectral space.

  • L (int, optional) – Denotes the period of the function to be approximated for the Fourier transform.

  • comm (COMM_WORLD, optional) – Communicator for mpi4py-fft.

fft

Object for parallel FFT transforms.

Type:

PFFT

X

Grid coordinates in real space.

Type:

mesh-grid

ndim

Number of spatial dimensions.

Type:

int

Ku

Laplace operator in spectral space (u component).

Type:

matrix

Kv

Laplace operator in spectral space (v component).

Type:

matrix

References

eval_f(u, t)[source]

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

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

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

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

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

Simple FFT solver for the diffusion part.

Parameters:
  • rhs (dtype_f) – Right-hand side for the linear 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 (e.g. for time-dependent BCs).

Returns:

me – Solution.

Return type:

dtype_u

u_exact(t)[source]

Routine to compute the exact solution at time \(t = 0\), see [3].

Parameters:

t (float) – Time of the exact solution.

Returns:

me – Exact solution.

Return type:

dtype_u

class grayscott_imex_linear(**kwargs)[source]

Bases: grayscott_imex_diffusion

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\). The model with linear (reaction) part is described by the \(N\)-dimensional model

\[\frac{d u}{d t} = D_u \Delta u - u v^2 + A,\]
\[\frac{d v}{d t} = D_v \Delta v + u v^2\]

in \(x \in \Omega:=[-L/2, L/2]^N\) with \(N=2,3\). Spatial discretization is done by using Fast Fourier transformation for solving the linear parts provided by mpi4py-fft [2]_, see also https://mpi4py-fft.readthedocs.io/en/latest/.

This class implements the problem for semi-explicit time-stepping (diffusion is treated implicitly, and linear part is computed in an explicit way).

eval_f(u, t)[source]

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

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

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

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

class grayscott_mi_diffusion(newton_maxiter=100, newton_tol=1e-12, **kwargs)[source]

Bases: grayscott_imex_diffusion

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\). Here, the process is described by the \(N\)-dimensional model

\[\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\]

in \(x \in \Omega:=[-L/2, L/2]^N\) with \(N=2,3\). Spatial discretization is done by using Fast Fourier transformation for solving the linear parts provided by mpi4py-fft [2]_, see also https://mpi4py-fft.readthedocs.io/en/latest/.

This class implements the problem for multi-implicit time-stepping, i.e., both diffusion and reaction part will be treated implicitly.

Parameters:
  • nvars (tuple of int, optional) – Spatial resolution, i.e., number of degrees of freedom in space. Should be a tuple, e.g. nvars=(127, 127).

  • 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\).

  • spectral (bool, optional) – If True, the solution is computed in spectral space.

  • L (int, optional) – Denotes the period of the function to be approximated for the Fourier transform.

  • comm (COMM_WORLD, optional) – Communicator for mpi4py-fft.

fft

Object for parallel FFT transforms.

Type:

PFFT

X

Grid coordinates in real space.

Type:

mesh-grid

ndim

Number of spatial dimensions.

Type:

int

Ku

Laplace operator in spectral space (u component).

Type:

matrix

Kv

Laplace operator in spectral space (v component).

Type:

matrix

References

dtype_f

alias of comp2_mesh

eval_f(u, t)[source]

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

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

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

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

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

Solver for the first component, can just call the super function.

Parameters:
  • rhs (dtype_f) – Right-hand side for the linear 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 (e.g. for time-dependent BCs).

Returns:

me – The solution as mesh.

Return type:

dtype_u

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

Newton-Solver for the second component.

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

  • float (factor) – 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 (e.g. for time-dependent BCs).

Returns:

me – The solution as mesh.

Return type:

dtype_u

class grayscott_mi_linear(newton_maxiter=100, newton_tol=1e-12, **kwargs)[source]

Bases: grayscott_imex_linear

The original 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\). The model with linear (reaction) part is described by the \(N\)-dimensional model

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

in \(x \in \Omega:=[-L/2, L/2]^N\) with \(N=2,3\). Spatial discretization is done by using Fast Fourier transformation for solving the linear parts provided by mpi4py-fft [2]_, see also https://mpi4py-fft.readthedocs.io/en/latest/.

The problem in this class will be treated in a multi-implicit way for time-stepping, i.e., for the system containing the diffusion part will be solved by FFT, and for the linear part a Newton solver is used.

dtype_f

alias of comp2_mesh

eval_f(u, t)[source]

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

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

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

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

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

Solver for the first component, can just call the super function.

Parameters:
  • rhs (dtype_f) – Right-hand side for the linear 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 (e.g. for time-dependent BCs).

Returns:

me – The solution as mesh.

Return type:

dtype_u

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

Newton-Solver for the second component.

Parameters:
  • rhs (dtype_f) – Right-hand side for the linear 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 (e.g. for time-dependent BCs).

Returns:

u – The solution as mesh.

Return type:

dtype_u