implementations.problem_classes.GrayScott_2D_PETSc_periodic module

class GS_full(da, prob, factor, dx, dy)[source]

Bases: object

Helper class to generate residual and Jacobian matrix for PETSc’s nonlinear solver SNES.

Parameters:
  • da (DMDA object) – Object of PETSc.

  • prob (problem instance) – Contains problem information for PETSc.

  • factor (float) – Temporal factor \(\Delta t Q_\Delta\).

  • dx (float) – Grid spacing in x direction.

  • dy (float) – Grid spacing in y direction.

localX

Local vector for PETSc.

Type:

PETSc vector object

formFunction(snes, X, F)[source]

Function to evaluate the residual for the Newton solver. This function should be equal to the RHS in the solution.

Parameters:
  • snes (PETSc solver object) – Nonlinear solver object.

  • X (PETSc vector object) – Input vector.

  • F (PETSc vector object) – Output vector \(F(X)\).

Returns:

Overwrites F.

Return type:

None

formJacobian(snes, X, J, P)[source]

Function to return the Jacobian matrix.

Parameters:
  • snes (PETSc solver object) – Nonlinear solver object.

  • X (PETSc vector object) – Input vector.

  • J (PETSc matrix object) – Current Jacobian matrix.

  • P (PETSc matrix object) – New Jacobian matrix.

Returns:

Matrix status.

Return type:

PETSc.Mat.Structure.SAME_NONZERO_PATTERN

class GS_reaction(da, prob, factor)[source]

Bases: object

Helper class to generate residual and Jacobian matrix for PETSc’s nonlinear solver SNES.

Parameters:
  • da (DMDA object) – Object of PETSc.

  • prob (problem instance) – Contains problem information for PETSc.

  • factor (float) – Temporal factor \(\Delta t Q_\Delta\).

localX

Local vector for PETSc.

Type:

PETSc vector object

formFunction(snes, X, F)[source]

Function to evaluate the residual for the Newton solver. This function should be equal to the RHS in the solution.

Parameters:
  • snes (PETSc solver object) – Nonlinear solver object.

  • X (PETSc vector object) – Input vector.

  • F (PETSc vector object) – Output vector \(F(X)\).

Returns:

Overwrites F.

Return type:

None

formJacobian(snes, X, J, P)[source]

Function to return the Jacobian matrix.

Parameters:
  • snes (PETSc solver object) – Nonlinear solver object.

  • X (PETSc vector object) – Input vector.

  • J (PETSc matrix object) – Current Jacobian matrix.

  • P (PETSc matrix object) – New Jacobian matrix.

Returns:

Matrix status.

Return type:

PETSc.Mat.Structure.SAME_NONZERO_PATTERN

class petsc_grayscott_fullyimplicit(nvars, Du, Dv, A, B, comm=petsc4py.PETSc.COMM_WORLD, lsol_tol=1e-10, nlsol_tol=1e-10, lsol_maxiter=None, nlsol_maxiter=None)[source]

Bases: petsc_grayscott_multiimplicit

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 two-dimensional model using periodic 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 spatial solve of the problem is realized by PETSc [2], [3]. For time-stepping, the problem is handled in a fully-implicit way, i.e., the nonlinear system containing the full right-hand side will be solved by PETSc’s nonlinear solver.

dtype_f

alias of petsc_vec

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 the numerical solution is computed.

Returns:

f – Right-hand side of the problem.

Return type:

dtype_f

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

Nonlinear solver for \((I - factor \cdot F)(\vec{u}) = \vec{rhs}\).

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

  • factor (float) – Abbrev. for the local stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver.

  • t (float) – Current time (e.g. for time-dependent BCs).

Returns:

me – Solution as mesh.

Return type:

dtype_u

class petsc_grayscott_multiimplicit(nvars, Du, Dv, A, B, comm=petsc4py.PETSc.COMM_WORLD, lsol_tol=1e-10, nlsol_tol=1e-10, lsol_maxiter=None, nlsol_maxiter=None)[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 two-dimensional model using periodic 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 spatial solve of the problem is realized by PETSc [2], [3]. For time-stepping, the diffusion part is solved by one of PETSc’s linear solver, whereas the reaction part will be solved by a nonlinear solver.

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

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

  • comm (PETSc.COMM_WORLD, optional) – Communicator for PETSc.

  • lsol_tol (float, optional) – Tolerance for linear solver to terminate.

  • nlsol_tol (float, optional) – Tolerance for nonlinear solver to terminate.

  • lsol_maxiter (int, optional) – Maximum number of iterations for linear solver.

  • nlsol_maxiter (int, optional) – Maximum number of iterations for nonlinear solver.

dx

Mesh grid width in x direction.

Type:

float

dy

Mesh grid width in y direction.

Type:

float

AMat

Discretization matrix.

Type:

PETSc matrix object

Id

Identity matrix.

Type:

PETSc matrix object

localX

Local vector for solution.

Type:

PETSc vector object

ksp

Linear solver object.

Type:

PETSc solver object

snes

Nonlinear solver object.

Type:

PETSc solver object

snes_itercount

Number of iterations done by nonlinear solver.

Type:

int

snes_ncalls

Number of calls of PETSc’s nonlinear solver.

Type:

int

References

dtype_f

alias of petsc_vec_comp2

dtype_u

alias of petsc_vec

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 the numerical solution is computed.

Returns:

f – Right-hand side of the problem.

Return type:

dtype_f

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

Linear solver for \((I - factor \cdot A)\vec{u} = \vec{rhs}\).

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

  • factor (float) – Abbrev. for the local stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver.

  • t (float) – Current time (e.g. for time-dependent BCs).

Returns:

me – Solution as mesh.

Return type:

dtype_u

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

Nonlinear solver for \((I - factor \cdot F)(\vec{u}) = \vec{rhs}\).

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

  • factor (float) – Abbrev. for the local stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver.

  • t (float) – Current time (e.g. for time-dependent BCs).

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.

Return type:

dtype_u

class petsc_grayscott_semiimplicit(nvars, Du, Dv, A, B, comm=petsc4py.PETSc.COMM_WORLD, lsol_tol=1e-10, nlsol_tol=1e-10, lsol_maxiter=None, nlsol_maxiter=None)[source]

Bases: petsc_grayscott_multiimplicit

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 two-dimensional model using periodic 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 spatial solve of the problem is realized by PETSc [2], [3]. For time-stepping, the problem is treated semi-implicitly, i.e., the system with diffusion part is solved by PETSc’s linear solver.

dtype_f

alias of petsc_vec_imex

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 the numerical solution is computed.

Returns:

f – Right-hand side of the problem.

Return type:

dtype_f

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

Linear solver for \((I - factor \cdot A)\vec{u} = \vec{rhs}\).

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

  • factor (float) – Abbrev. for the local stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver.

  • t (float) – Current time (e.g. for time-dependent BCs).

Returns:

me – Solution as mesh.

Return type:

dtype_u