HyPar
1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
|
Contains the functions required for Jacobian computations for implicit time-integration. More...
#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <mpivars_cpp.h>
#include <simulation_object.h>
#include <petscinterface.h>
Go to the source code of this file.
Macros | |
#define | __FUNCT__ "PetscIJacobian" |
#define | __FUNCT__ "PetscJacobianFunction_JFNK" |
#define | __FUNCT__ "PetscJacobianFunction_Linear" |
Functions | |
PetscErrorCode | PetscIJacobian (TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *ctxt) |
PetscErrorCode | PetscJacobianFunction_JFNK (Mat Jacobian, Vec Y, Vec F) |
PetscErrorCode | PetscJacobianFunction_Linear (Mat Jacobian, Vec Y, Vec F) |
Contains the functions required for Jacobian computations for implicit time-integration.
Definition in file PetscIJacobian.cpp.
#define __FUNCT__ "PetscIJacobian" |
Definition at line 188 of file PetscIJacobian.cpp.
#define __FUNCT__ "PetscJacobianFunction_JFNK" |
Definition at line 188 of file PetscIJacobian.cpp.
#define __FUNCT__ "PetscJacobianFunction_Linear" |
Definition at line 188 of file PetscIJacobian.cpp.
PetscErrorCode PetscIJacobian | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | Ydot, | ||
PetscReal | a, | ||
Mat | A, | ||
Mat | B, | ||
void * | ctxt | ||
) |
Compute the Jacobian for implicit time integration of the governing equations: The ODE, obtained after discretizing the governing PDE in space, is expressed as follows:
\begin{equation} \frac {d{\bf U}}{dt} = {\bf F}\left({\bf U}\right) \Rightarrow \dot{\bf U} - {\bf F}\left({\bf U}\right) = 0, \end{equation}
where {\bf F} is the spatially discretized right-hand-side, and {\bf U} represents the entire solution vector. The Jacobian is thus given by:
\begin{equation} {\bf J} = \left[\alpha{\bf I} - \frac {\partial {\bf F}} {\partial {\bf U}} \right] \end{equation}
where \alpha is the shift coefficient (PETScContext::shift) of the time integration method.
Matrix-free representation: The Jacobian is computed using a matrix-free approach, where the entire Jacobian is never assembled, computed, or stored. Instead, the action of the Jacobian on a vector is defined through functions (PetscJacobianFunction_JFNK() for nonlinear systems, and PetscJacobianFunction_Linear() for linear systems). This approach works well with iterative solvers. Thus, this function just does the following:
Notes:
ts | Time stepping object (see PETSc TS) |
t | Current time |
Y | Solution vector |
Ydot | Time-derivative of solution vector |
a | Shift |
A | Jacobian matrix |
B | Preconditioning matrix |
ctxt | Application context |
Definition at line 52 of file PetscIJacobian.cpp.
PetscErrorCode PetscJacobianFunction_JFNK | ( | Mat | Jacobian, |
Vec | Y, | ||
Vec | F | ||
) |
Computes the action of the Jacobian on a vector: See documentation for PetscIJacobian() for the definition of the Jacobian. This function computes its action on a vector for a nonlinear system by taking the directional derivative, as follows:
\begin{equation} {\bf f} = {\bf J}\left({\bf U}_0\right){\bf y} = \frac {\partial \mathcal{F}\left({\bf U}\right)} {\partial {\bf U}}\left({\bf U}_0\right) {\bf y} \approx \frac{1}{\epsilon} \left[ \mathcal{F}\left({\bf U}_0+\epsilon{\bf y}\right)-\mathcal{F}\left({\bf U}_0\right) \right] \end{equation}
In the context of the implicit time integration of the ODE given by
\begin{equation} \frac {d {\bf U}} {dt} = {\bf F}\left({\bf U}\right) \Rightarrow \frac {d {\bf U}} {dt} - {\bf F}\left({\bf U}\right) = 0, \end{equation}
we have that
\begin{equation} \mathcal{F}\left({\bf U}\right) \equiv \dot{\bf U} - {\bf F}\left({\bf U}\right) \Rightarrow {\bf J}\left({\bf U}\right) = \left[\alpha{\bf I} - \frac {\partial {\bf F}\left({\bf U}\right)} {\partial {\bf U}} \right]. \end{equation}
where \alpha (PETScContext::shift) is the shift variable specific to the time integration method. So this function computes
\begin{equation} {\bf f} = \alpha {\bf y} - \frac{1}{\epsilon} \left[ {\bf F}\left({\bf U}_0+\epsilon{\bf y}\right)-{\bf F}\left({\bf U}_0\right) \right] \end{equation}
In the code, {\bf y} is Y, \bf f is F, and {\bf U}_0 is HyPar::uref (the reference solution at which the nonlinear Jacobian is computed). See papers on Jacobian-free Newton-Krylov (JFNK) methods to understand how \epsilon is computed.
Notes:
Jacobian | Jacobian matrix |
Y | Input vector |
F | Output vector (Jacobian times input vector) |
Definition at line 111 of file PetscIJacobian.cpp.
PetscErrorCode PetscJacobianFunction_Linear | ( | Mat | Jacobian, |
Vec | Y, | ||
Vec | F | ||
) |
Computes the action of the Jacobian on a vector: See documentation for PetscIJacobian() for the definition of the Jacobian. This function computes its action on a vector for a linear system by taking the directional derivative, as follows:
\begin{equation} {\bf f} = {\bf J}{\bf y} = \frac {\partial \mathcal{F}\left({\bf U}\right)} {\partial {\bf U}} {\bf y} = \left[ \mathcal{F}\left({\bf U}_0+{\bf y}\right)-\mathcal{F}\left({\bf U}_0\right) \right], \end{equation}
where \mathcal{F} is linear, and thus {\bf J} is a constant ( \mathcal{F}\left({\bf y}\right) = {\bf J}{\bf y}). In the context of the implicit time integration of a linear ODE given by
\begin{equation} \frac {d {\bf U}} {dt} = {\bf F}\left({\bf U}\right) \Rightarrow \frac {d {\bf U}} {dt} - {\bf F}\left({\bf U}\right) = 0, \end{equation}
we have that
\begin{equation} \mathcal{F}\left({\bf U}\right) \equiv \dot{\bf U} - {\bf F}\left({\bf U}\right) \Rightarrow {\bf J}\left({\bf U}\right) = \left[\alpha{\bf I} - \frac {\partial {\bf F}\left({\bf U}\right)} {\partial {\bf U}} \right]. \end{equation}
where \alpha (PETScContext::shift) is the shift variable specific to the time integration method. So this function computes
\begin{equation} {\bf f} = \alpha {\bf y} - \left[ {\bf F}\left({\bf U}_0+{\bf y}\right)-{\bf F}\left({\bf U}_0\right) \right] \end{equation}
In the code, {\bf y} is Y, \bf f is F, and {\bf U}_0 is HyPar::uref (the reference solution at which the nonlinear Jacobian is computed).
Since \mathcal{F} is linear,
\begin{equation} {\bf J}{\bf y} = \left[ \mathcal{F}\left({\bf U}_0+{\bf y}\right)-\mathcal{F}\left({\bf U}_0\right) \right] = \mathcal{F}\left({\bf y}\right). \end{equation}
However, the Jacobian is not computed as \mathcal{F}\left({\bf y}\right) because of the following reason: this function is used by the equation solver within the implicit time integrator in PETSc, and {\bf y} represents the change in the solution, i.e. \Delta {\bf U}, and not the solution \bf U. Thus, evaluating \mathcal{F}\left({\bf y}\right) using HyPar::HyperbolicFunction, HyPar::ParabolicFunction, and HyPar::SourceFunction is incorrect since these functions expect \bf U as the input.
Notes:
Jacobian | Jacobian matrix |
Y | Input vector |
F | Output vector (Jacobian times input vector |
Definition at line 236 of file PetscIJacobian.cpp.