HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
PetscIJacobian.cpp File Reference

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)
 

Detailed Description

Contains the functions required for Jacobian computations for implicit time-integration.

Author
Debojyoti Ghosh

Definition in file PetscIJacobian.cpp.

Macro Definition Documentation

◆ __FUNCT__ [1/3]

#define __FUNCT__   "PetscIJacobian"

Definition at line 188 of file PetscIJacobian.cpp.

◆ __FUNCT__ [2/3]

#define __FUNCT__   "PetscJacobianFunction_JFNK"

Definition at line 188 of file PetscIJacobian.cpp.

◆ __FUNCT__ [3/3]

#define __FUNCT__   "PetscJacobianFunction_Linear"

Definition at line 188 of file PetscIJacobian.cpp.

Function Documentation

◆ PetscIJacobian()

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:

Parameters
tsTime stepping object (see PETSc TS)
tCurrent time
YSolution vector
YdotTime-derivative of solution vector
aShift
AJacobian matrix
BPreconditioning matrix
ctxtApplication context

Definition at line 52 of file PetscIJacobian.cpp.

60 {
61  PETScContext* context = (PETScContext*) ctxt;
62 
63  PetscFunctionBegin;
64  for (int ns = 0; ns < context->nsims; ns++) {
65  ((SimulationObject*)context->simobj)[ns].solver.count_IJacobian++;
66  }
67  context->shift = a;
68  context->waqt = t;
69  /* Construct preconditioning matrix */
70  if (context->flag_use_precon) PetscComputePreconMatImpl(B,Y,context);
71 
72  PetscFunctionReturn(0);
73 }
Structure defining a simulation.
int PetscComputePreconMatImpl(Mat, Vec, void *)
Structure containing the variables for time-integration with PETSc.

◆ PetscJacobianFunction_JFNK()

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:

  • For a nonlinear spatial discretization, the right-hand-side must be computed without the nonlinearity (i.e. with a previously computed or "frozen" discretization operator). This ensures that the Jacobian being computed is consistent.
  • See https://petsc.org/release/docs/manualpages/TS/index.html for documentation on PETSc's time integrators.
  • All functions and variables whose names start with Vec, Mat, PC, KSP, SNES, and TS are defined by PETSc. Refer to the PETSc documentation (https://petsc.org/release/docs/). Usually, googling with the function or variable name yields the specific doc page dealing with that function/variable.
Parameters
JacobianJacobian matrix
YInput vector
FOutput vector (Jacobian times input vector)

Definition at line 111 of file PetscIJacobian.cpp.

115 {
116  PETScContext* context(nullptr);
117 
118  PetscFunctionBegin;
119 
120  MatShellGetContext(Jacobian,&context);
121  SimulationObject* sim = (SimulationObject*) context->simobj;
122  int nsims = context->nsims;
123 
124  double normY;
125  VecNorm(Y,NORM_2,&normY);
126 
127  if (normY < 1e-16) {
128 
129  /* F = 0 */
130  VecZeroEntries(F);
131  /* [J]Y = aY - F(Y) */
132  VecAXPBY(F,context->shift,0,Y);
133 
134  } else {
135 
136  double epsilon = context->jfnk_eps / normY;
137  double t = context->waqt; /* current stage/step time */
138 
139  for (int ns = 0; ns < nsims; ns++) {
140 
141  HyPar* solver = &(sim[ns].solver);
142  MPIVariables* mpi = &(sim[ns].mpi);
143  solver->count_IJacFunction++;
144 
145  int size = solver->npoints_local_wghosts;
146 
147  double *u = solver->u;
148  double *uref = solver->uref;
149  double *rhsref = solver->rhsref;
150  double *rhs = solver->rhs;
151 
152  /* copy solution from PETSc vector */
153  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
154  _ArrayAYPX_(uref,epsilon,u,size*solver->nvars);
155  /* apply boundary conditions and exchange data over MPI interfaces */
156  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
158  solver->nvars,
159  solver->dim_local,
160  solver->ghosts,
161  mpi,
162  u );
163 
164  /* Evaluate hyperbolic, parabolic and source terms and the RHS for U+dU */
165  _ArraySetValue_(rhs,size*solver->nvars,0.0);
166  solver->HyperbolicFunction( solver->hyp,u,solver,mpi,t,0,
167  solver->FFunction,solver->Upwind);
168  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
169  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
170  _ArrayAXPY_(solver->par, 1.0,rhs,size*solver->nvars);
171  solver->SourceFunction (solver->source,u,solver,mpi,t);
172  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
173 
174  _ArrayAXPY_(rhsref,-1.0,rhs,size*solver->nvars);
175  /* Transfer RHS to PETSc vector */
176  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
177  }
178 
179  /* [J]Y = aY - F(Y) */
180  VecAXPBY(F,context->shift,(-1.0/epsilon),Y);
181 
182  }
183 
184  PetscFunctionReturn(0);
185 }
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
int nvars
Definition: hypar.h:29
int count_IJacFunction
Definition: hypar.h:422
int(* ParabolicFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:256
Structure defining a simulation.
double * par
Definition: hypar.h:122
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
int(* ApplyBoundaryConditions)(void *, void *, double *, double *, double)
Definition: hypar.h:214
double * u
Definition: hypar.h:116
#define _ArrayAYPX_(x, a, y, size)
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int(* SourceFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:259
int ndims
Definition: hypar.h:26
int(* HyperbolicFunction)(double *, double *, void *, void *, double, int, int(*)(double *, double *, int, void *, double), int(*)(double *, double *, double *, double *, double *, double *, int, void *, double))
Definition: hypar.h:250
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
double * rhs
Definition: hypar.h:399
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArrayAXPY_(x, a, y, size)
Structure containing the variables for time-integration with PETSc.
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
double * rhsref
Definition: hypar.h:398
int ghosts
Definition: hypar.h:52
Structure of MPI-related variables.
int npoints_local_wghosts
Definition: hypar.h:42
double * uref
Definition: hypar.h:397
int TransferVecFromPETSc(double *const, const Vec, void *, const int, const int)
double * source
Definition: hypar.h:125
double * hyp
Definition: hypar.h:119

◆ PetscJacobianFunction_Linear()

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:

  • For a nonlinear spatial discretization, the right-hand-side must be computed without the nonlinearity (i.e. with a previously computed or "frozen" discretization operator). This ensures that the Jacobian being computed is consistent, and is truly linear.
  • See https://petsc.org/release/docs/manualpages/TS/index.html for documentation on PETSc's time integrators.
  • All functions and variables whose names start with Vec, Mat, PC, KSP, SNES, and TS are defined by PETSc. Refer to the PETSc documentation (https://petsc.org/release/docs/). Usually, googling with the function or variable name yields the specific doc page dealing with that function/variable.
Parameters
JacobianJacobian matrix
YInput vector
FOutput vector (Jacobian times input vector

Definition at line 236 of file PetscIJacobian.cpp.

240 {
241  PETScContext* context(nullptr);
242 
243  PetscFunctionBegin;
244 
245  MatShellGetContext(Jacobian,&context);
246  SimulationObject* sim = (SimulationObject*) context->simobj;
247  int nsims = context->nsims;
248 
249  double normY;
250  VecNorm(Y,NORM_2,&normY);
251 
252  if (normY < 1e-16) {
253 
254  /* F = 0 */
255  VecZeroEntries(F);
256  /* [J]Y = aY - F(Y) */
257  VecAXPBY(F,context->shift,0,Y);
258 
259  } else {
260 
261  double t = context->waqt; /* current stage/step time */
262 
263  for (int ns = 0; ns < nsims; ns++) {
264 
265  HyPar* solver = &(sim[ns].solver);
266  MPIVariables* mpi = &(sim[ns].mpi);
267  solver->count_IJacFunction++;
268 
269  int size = solver->npoints_local_wghosts;
270 
271  double *u = solver->u;
272  double *uref = solver->uref;
273  double *rhsref = solver->rhsref;
274  double *rhs = solver->rhs;
275 
276  /* copy solution from PETSc vector */
277  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
278  _ArrayAYPX_(uref,1.0,u,size*solver->nvars);
279  /* apply boundary conditions and exchange data over MPI interfaces */
280  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
282  solver->nvars,
283  solver->dim_local,
284  solver->ghosts,
285  mpi,
286  u );
287 
288  /* Evaluate hyperbolic, parabolic and source terms and the RHS for U+dU */
289  _ArraySetValue_(rhs,size*solver->nvars,0.0);
290  solver->HyperbolicFunction( solver->hyp,u,solver,mpi,t,0,
291  solver->FFunction,solver->Upwind);
292  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
293  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
294  _ArrayAXPY_(solver->par, 1.0,rhs,size*solver->nvars);
295  solver->SourceFunction (solver->source,u,solver,mpi,t);
296  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
297 
298  _ArrayAXPY_(rhsref,-1.0,rhs,size*solver->nvars);
299  /* Transfer RHS to PETSc vector */
300  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
301  }
302 
303  /* [J]Y = aY - F(Y) */
304  VecAXPBY(F,context->shift,-1.0,Y);
305 
306  }
307 
308  PetscFunctionReturn(0);
309 }
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
int nvars
Definition: hypar.h:29
int count_IJacFunction
Definition: hypar.h:422
int(* ParabolicFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:256
Structure defining a simulation.
double * par
Definition: hypar.h:122
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
int(* ApplyBoundaryConditions)(void *, void *, double *, double *, double)
Definition: hypar.h:214
double * u
Definition: hypar.h:116
#define _ArrayAYPX_(x, a, y, size)
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int(* SourceFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:259
int ndims
Definition: hypar.h:26
int(* HyperbolicFunction)(double *, double *, void *, void *, double, int, int(*)(double *, double *, int, void *, double), int(*)(double *, double *, double *, double *, double *, double *, int, void *, double))
Definition: hypar.h:250
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
double * rhs
Definition: hypar.h:399
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArrayAXPY_(x, a, y, size)
Structure containing the variables for time-integration with PETSc.
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
double * rhsref
Definition: hypar.h:398
int ghosts
Definition: hypar.h:52
Structure of MPI-related variables.
int npoints_local_wghosts
Definition: hypar.h:42
double * uref
Definition: hypar.h:397
int TransferVecFromPETSc(double *const, const Vec, void *, const int, const int)
double * source
Definition: hypar.h:125
double * hyp
Definition: hypar.h:119