implementations.problem_classes.AllenCahn_2D_FD module¶
- class allencahn_fullyimplicit(nvars=(128, 128), nu=2, eps=0.04, newton_maxiter=200, newton_tol=1e-12, lin_tol=1e-08, lin_maxiter=100, inexact_linear_ratio=None, radius=0.25, order=2)[source]¶
Bases:
ptype
Example implementing the two-dimensional Allen-Cahn equation with periodic boundary conditions \(u \in [-1, 1]^2\)
\[\frac{\partial u}{\partial t} = \Delta u + \frac{1}{\varepsilon^2} u (1 - u^\nu)\]for constant parameter \(\nu\). Initial condition are circles of the form
\[u({\bf x}, 0) = \tanh\left(\frac{r - \sqrt{x_i^2 + y_j^2}}{\sqrt{2}\varepsilon}\right)\]for \(i, j=0,..,N-1\), where \(N\) is the number of spatial grid points. For time-stepping, the problem is treated fully-implicitly, i.e., the nonlinear system is solved by Newton.
- Parameters:
nvars (tuple of int, optional) – Number of unknowns in the problem, e.g.
nvars=(128, 128)
.nu (float, optional) – Problem parameter \(\nu\).
eps (float, optional) – Scaling parameter \(\varepsilon\).
newton_maxiter (int, optional) – Maximum number of iterations for the Newton solver.
newton_tol (float, optional) – Tolerance for Newton’s method to terminate.
lin_tol (float, optional) – Tolerance for linear solver to terminate.
lin_maxiter (int, optional) – Maximum number of iterations for the linear solver.
radius (float, optional) – Radius of the circles.
order (int, optional) – Order of the finite difference matrix.
- A¶
Second-order FD discretization of the 2D laplace operator.
- Type:
scipy.spdiags
- dx¶
Distance between two spatial nodes (same for both directions).
- Type:
float
- xvalues¶
Spatial grid points, here both dimensions have the same grid points.
- Type:
np.1darray
- newton_itercount¶
Number of iterations of Newton solver.
- Type:
int
- lin_itercount¶
Number of iterations of linear solver.
- newton_ncalls¶
Number of calls of Newton solver.
- Type:
int
- lin_ncalls¶
Number of calls of linear solver.
- Type:
int
- dtype_f¶
alias of
mesh
- dtype_u¶
alias of
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 (not used here).
- Returns:
f – The right-hand side of the problem.
- Return type:
dtype_f
- solve_system(rhs, factor, u0, t)[source]¶
Simple Newton solver.
- 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.
t (float) – Current time (required here for the BC).
- Returns:
me – The solution as mesh.
- Return type:
dtype_u
- class allencahn_multiimplicit(nvars=(128, 128), nu=2, eps=0.04, newton_maxiter=200, newton_tol=1e-12, lin_tol=1e-08, lin_maxiter=100, inexact_linear_ratio=None, radius=0.25, order=2)[source]¶
Bases:
allencahn_fullyimplicit
Example implementing the two-dimensional Allen-Cahn equation with periodic boundary conditions \(u \in [-1, 1]^2\)
\[\frac{\partial u}{\partial t} = \Delta u + \frac{1}{\varepsilon^2} u (1 - u^\nu)\]for constant parameter \(\nu\). Initial condition are circles of the form
\[u({\bf x}, 0) = \tanh\left(\frac{r - \sqrt{x_i^2 + y_j^2}}{\sqrt{2}\varepsilon}\right)\]for \(i, j=0,..,N-1\), where \(N\) is the number of spatial grid points. For time-stepping, the problem is treated in multi-implicit fashion, i.e., the linear system containing the Laplacian is solved by the conjugate gradients method, and the system containing the rest of the right-hand side will be solved by Newton’s method.
- 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]¶
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 – The solution as mesh.
- Return type:
dtype_u
- solve_system_2(rhs, factor, u0, t)[source]¶
Simple Newton solver.
- 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.
t (float) – Current time (required here for the BC).
- Returns:
me – The solution as mesh.
- Return type:
dtype_u
- class allencahn_multiimplicit_v2(nvars=(128, 128), nu=2, eps=0.04, newton_maxiter=200, newton_tol=1e-12, lin_tol=1e-08, lin_maxiter=100, inexact_linear_ratio=None, radius=0.25, order=2)[source]¶
Bases:
allencahn_fullyimplicit
This class implements the two-dimensional Allen-Cahn (AC) equation with periodic boundary conditions \(u \in [-1, 1]^2\)
\[\frac{\partial u}{\partial t} = \Delta u + \frac{1}{\varepsilon^2} u (1 - u^\nu)\]for constant parameter \(\nu\). The initial condition has the form of circles
\[u({\bf x}, 0) = \tanh\left(\frac{r - \sqrt{x_i^2 + y_j^2}}{\sqrt{2}\varepsilon}\right)\]for \(i, j=0,..,N-1\), where \(N\) is the number of spatial grid points. For time-stepping, a special AC-splitting is used here to get another kind of semi-implicit treatment of the problem: The term \(\Delta u - \frac{1}{\varepsilon^2} u^{\nu + 1}\) is handled implicitly and the nonlinear system including this part will be solved by Newton. \(\frac{1}{\varepsilon^2} u\) is solved by a linear solver provided by a
SciPy
routine.- 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]¶
Simple Newton solver.
- 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.
t (float) – Current time (required here for the BC).
- Returns:
me – The solution as mesh.
- Return type:
dtype_u
- solve_system_2(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 – The solution as mesh.
- Return type:
dtype_u
- class allencahn_semiimplicit(nvars=(128, 128), nu=2, eps=0.04, newton_maxiter=200, newton_tol=1e-12, lin_tol=1e-08, lin_maxiter=100, inexact_linear_ratio=None, radius=0.25, order=2)[source]¶
Bases:
allencahn_fullyimplicit
This class implements the two-dimensional Allen-Cahn equation with periodic boundary conditions \(u \in [-1, 1]^2\)
\[\frac{\partial u}{\partial t} = \Delta u + \frac{1}{\varepsilon^2} u (1 - u^\nu)\]for constant parameter \(\nu\). Initial condition are circles of the form
\[u({\bf x}, 0) = \tanh\left(\frac{r - \sqrt{x_i^2 + y_j^2}}{\sqrt{2}\varepsilon}\right)\]for \(i, j=0,..,N-1\), where \(N\) is the number of spatial grid points. For time-stepping, the problem is treated in a semi-implicit way, i.e., the linear system containing the Laplacian is solved by the conjugate gradients method, and the system containing the rest of the right-hand side is only evaluated at each time.
- dtype_f¶
alias of
imex_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 (not used here).
- Returns:
f – The 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 – The solution as mesh.
- Return type:
dtype_u
- class allencahn_semiimplicit_v2(nvars=(128, 128), nu=2, eps=0.04, newton_maxiter=200, newton_tol=1e-12, lin_tol=1e-08, lin_maxiter=100, inexact_linear_ratio=None, radius=0.25, order=2)[source]¶
Bases:
allencahn_fullyimplicit
This class implements the two-dimensional Allen-Cahn (AC) equation with periodic boundary conditions \(u \in [-1, 1]^2\)
\[\frac{\partial u}{\partial t} = \Delta u + \frac{1}{\varepsilon^2} u (1 - u^\nu)\]for constant parameter \(\nu\). Initial condition are circles of the form
\[u({\bf x}, 0) = \tanh\left(\frac{r - \sqrt{x_i^2 + y_j^2}}{\sqrt{2}\varepsilon}\right)\]for \(i, j=0,..,N-1\), where \(N\) is the number of spatial grid points. For time-stepping, a special AC-splitting is used to get a semi-implicit treatment of the problem: The term \(\Delta u - \frac{1}{\varepsilon^2} u^{\nu + 1}\) is handled implicitly and the nonlinear system including this part will be solved by Newton. \(\frac{1}{\varepsilon^2} u\) is only evaluated at each time.
- dtype_f¶
alias of
imex_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(rhs, factor, u0, t)[source]¶
Simple Newton solver.
- 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.
t (float) – Current time (required here for the BC).
- Returns:
me – The solution as mesh.
- Return type:
dtype_u