implementations.problem_classes.GeneralizedFisher_1D_PETSc module¶
- class Fisher_full(da, prob, factor, dx)[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.
- localX¶
Local vector for
PETSc
.- Type:
PETSc vector object
- xs, xe
Defines the ranges for spatial domain.
- Type:
int
- mx¶
Get sizes for the vector containing the spatial points.
- Type:
tuple
- row¶
Row for a matrix.
- Type:
PETSc matrix stencil object
- col¶
Column for a matrix.
- Type:
PETSc matrix stencil 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 Fisher_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\).
dx (float) – Grid spacing in x 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 petsc_fisher_fullyimplicit(nvars, lambda0, nu, interval, comm=petsc4py.PETSc.COMM_WORLD, lsol_tol=1e-10, nlsol_tol=1e-10, lsol_maxiter=None, nlsol_maxiter=None)[source]¶
Bases:
petsc_fisher_multiimplicit
The following one-dimensional problem is an example of a reaction-diffusion equation with traveling waves, and can be seen as a generalized Fisher equation. This class implements a special case of the Kolmogorov-Petrovskii-Piskunov problem [1] using periodic boundary conditions
\[\frac{\partial u}{\partial t} = \Delta u + \lambda_0^2 u (1 - u^\nu)\]with exact solution
\[u(x, 0) = \left[ 1 + \left(2^{\nu / 2} - 1\right) \exp\left(-(\nu / 2)\sigma_1 x + 2 \lambda_1 t\right) \right]^{-2 / \nu}\]for \(x \in \mathbb{R}\), and
\[\sigma_1 = \lambda_1 - \sqrt{\lambda_1^2 - \lambda_0^2},\]\[\lambda_1 = \frac{\lambda_0}{2} \left[ \left(1 + \frac{\nu}{2}\right)^{1/2} + \left(1 + \frac{\nu}{2}\right)^{-1/2} \right].\]This class is implemented to be solved in spatial using
PETSc
[2], [3]. For time-stepping, the problem is treated fully-implicitly.- 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_fisher_multiimplicit(nvars, lambda0, nu, interval, comm=petsc4py.PETSc.COMM_WORLD, lsol_tol=1e-10, nlsol_tol=1e-10, lsol_maxiter=None, nlsol_maxiter=None)[source]¶
Bases:
ptype
The following one-dimensional problem is an example of a reaction-diffusion equation with traveling waves, and can be seen as a generalized Fisher equation. This class implements a special case of the Kolmogorov-Petrovskii-Piskunov problem [1] using periodic boundary conditions
\[\frac{\partial u}{\partial t} = \Delta u + \lambda_0^2 u (1 - u^\nu)\]with exact solution
\[u(x, 0) = \left[ 1 + \left(2^{\nu / 2} - 1\right) \exp\left(-(\nu / 2)\sigma_1 x + 2 \lambda_1 t\right) \right]^{-2 / \nu}\]for \(x \in \mathbb{R}\), and
\[\sigma_1 = \lambda_1 - \sqrt{\lambda_1^2 - \lambda_0^2},\]\[\lambda_1 = \frac{\lambda_0}{2} \left[ \left(1 + \frac{\nu}{2}\right)^{1/2} + \left(1 + \frac{\nu}{2}\right)^{-1/2} \right].\]This class is implemented to be solved in spatial using
PETSc
[2], [3]. For time-stepping, the problem will be solved multi-implicitly.- Parameters:
nvars (int, optional) – Spatial resolution.
lambda0 (float, optional) – Problem parameter : math:lambda_0.
nu (float, optional) – Problem parameter \(\nu\).
interval (tuple, optional) – Defines the spatial domain.
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.
- Type:
float
- xs, xe
Define the ranges.
- Type:
int
- A¶
Discretization 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
- F¶
Global vector.
- Type:
PETSc vector object
- J¶
Jacobi matrix.
- Type:
PETSc matrix object
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
- get_sys_mat(factor)[source]¶
Helper function to assemble the system matrix of the linear problem.
- Parameters:
factor (float) – Factor to define the system matrix.
- Returns:
A – Matrix for the system to solve.
- Return type:
PETSc matrix object
- 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
- class petsc_fisher_semiimplicit(nvars, lambda0, nu, interval, comm=petsc4py.PETSc.COMM_WORLD, lsol_tol=1e-10, nlsol_tol=1e-10, lsol_maxiter=None, nlsol_maxiter=None)[source]¶
Bases:
petsc_fisher_multiimplicit
The following one-dimensional problem is an example of a reaction-diffusion equation with traveling waves, and can be seen as a generalized Fisher equation. This class implements a special case of the Kolmogorov-Petrovskii-Piskunov problem [1] using periodic boundary conditions
\[\frac{\partial u}{\partial t} = \Delta u + \lambda_0^2 u (1 - u^\nu)\]with exact solution
\[u(x, 0) = \left[ 1 + \left(2^{\nu / 2} - 1\right) \exp\left(-(\nu / 2)\sigma_1 x + 2 \lambda_1 t\right) \right]^{-2 / \nu}\]for \(x \in \mathbb{R}\), and
\[\sigma_1 = \lambda_1 - \sqrt{\lambda_1^2 - \lambda_0^2},\]\[\lambda_1 = \frac{\lambda_0}{2} \left[ \left(1 + \frac{\nu}{2}\right)^{1/2} + \left(1 + \frac{\nu}{2}\right)^{-1/2} \right].\]This class is implemented to be solved in spatial using
PETSc
[2], [3]. For time-stepping, the problem here will be solved in a semi-implicit way.- 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]¶
Simple 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