HyPar
1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
|
#include <petscinterface_struct.h>
Go to the source code of this file.
Functions | |
int | TransferVecToPETSc (const double *const, Vec, void *, const int, const int) |
int | TransferVecFromPETSc (double *const, const Vec, void *, const int, const int) |
int | TransferMatToPETSc (void *, Mat, void *) |
int | PetscRegisterTIMethods (int) |
PetscErrorCode | PetscRHSFunctionExpl (TS, PetscReal, Vec, Vec, void *) |
PetscErrorCode | PetscRHSFunctionIMEX (TS, PetscReal, Vec, Vec, void *) |
PetscErrorCode | PetscIFunctionIMEX (TS, PetscReal, Vec, Vec, Vec, void *) |
PetscErrorCode | PetscIFunctionImpl (TS, PetscReal, Vec, Vec, Vec, void *) |
PetscErrorCode | PetscIJacobianIMEX (TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *) |
PetscErrorCode | PetscJacobianFunctionIMEX_JFNK (Mat, Vec, Vec) |
PetscErrorCode | PetscJacobianFunctionIMEX_Linear (Mat, Vec, Vec) |
PetscErrorCode | PetscIJacobian (TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *) |
PetscErrorCode | PetscJacobianFunction_JFNK (Mat, Vec, Vec) |
PetscErrorCode | PetscJacobianFunction_Linear (Mat, Vec, Vec) |
int | PetscGlobalDOF (void *) |
int | PetscCleanup (void *) |
int | PetscCreatePointList (void *) |
PetscErrorCode | PetscPreStage (TS, PetscReal) |
PetscErrorCode | PetscPostStage (TS, PetscReal, PetscInt, Vec *) |
PetscErrorCode | PetscPreTimeStep (TS) |
PetscErrorCode | PetscPostTimeStep (TS) |
int | PetscComputePreconMatIMEX (Mat, Vec, void *) |
int | PetscComputePreconMatImpl (Mat, Vec, void *) |
int | PetscJacobianMatNonzeroEntriesImpl (Mat, int, void *) |
PetscErrorCode | PetscTimeError (TS) |
PetscErrorCode | PetscSetInitialGuessROM (SNES, Vec, void *) |
int TransferVecToPETSc | ( | const double *const | u, |
Vec | Y, | ||
void * | ctxt, | ||
const int | sim_idx, | ||
const int | offset | ||
) |
Copy data to a PETSc vector (used by PETSc time integrators, and with no ghost points) from a HyPar::u array (with ghost points).
u | HyPar::u type array (with ghost points) |
Y | PETSc vector |
ctxt | Object of type PETScContext |
sim_idx | Simulation object index |
offset | Offset |
Definition at line 24 of file TransferToPETSc.cpp.
int TransferVecFromPETSc | ( | double *const | u, |
const Vec | Y, | ||
void * | ctxt, | ||
const int | sim_idx, | ||
const int | offset | ||
) |
Copy data from a PETSc vector (used by PETSc time integrators, and with no ghost points) to a HyPar::u array (with ghost points).
u | HyPar::u type array (with ghost points) |
Y | PETSc vector |
ctxt | Object of type PETScContext |
sim_idx | Simulation object index |
offset | Offset |
Definition at line 23 of file TransferFromPETSc.cpp.
int TransferMatToPETSc | ( | void * | J, |
Mat | A, | ||
void * | ctxt | ||
) |
Copy a matrix of type BandedMatrix to a PETSc matrix.
J | Matrix of type BandedMatrix |
A | PETSc matrix |
ctxt | Object of type PETScContext |
Definition at line 53 of file TransferToPETSc.cpp.
int PetscRegisterTIMethods | ( | int | rank | ) |
This function allows the registering of custom multi-stage time integration methods and using it with PETSc. More than one method may be provided, each with a unique name (the name must not conflict with the names of existing methods in PETSc). The methods should be provided in a file "time_method.inp", and the following must be provided:
The method can then be used by invoking it by its name. For example, if a custom ARKIMEX method "foo" is provided through "time_method.inp", it can be used by specifying "-ts_type arkimex -ts_arkimex_type foo" in the command line or ".petscrc" file.
See /Examples/PETScInputs for examples of the input files required to provide a time integration method.
Currently supported classes of time integrators for which custom methods can be provided:
To do:
rank | MPI rank |
Definition at line 46 of file PetscRegisterTIMethods.cpp.
PetscErrorCode PetscRHSFunctionExpl | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | F, | ||
void * | ctxt | ||
) |
Compute the right-hand-side (RHS) for the explicit time integration of the governing equations: The spatially discretized ODE can be expressed as
\begin{equation} \frac {d{\bf U}} {dt} = {\bf F}\left({\bf U}\right). \end{equation}
This function computes \({\bf F}\left({\bf U}\right)\), given \({\bf U}\).
Note:
ts | Time integration object |
t | Current simulation time |
Y | State vector (input) |
F | The computed right-hand-side vector |
ctxt | Object of type PETScContext |
Definition at line 36 of file PetscRHSFunctionExpl.cpp.
PetscErrorCode PetscRHSFunctionIMEX | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | F, | ||
void * | ctxt | ||
) |
Compute the explicitly-treated right-hand-side (RHS) for the implicit-explicit (IMEX) time integration of the governing equations: The ODE, obtained after discretizing the governing PDE in space, is expressed as follows (for the purpose of IMEX time integration):
\begin{eqnarray} \frac {d{\bf U}}{dt} &=& {\bf F}\left({\bf U}\right) + {\bf G}\left({\bf U}\right), \\ \Rightarrow \dot{\bf U} - {\bf G}\left({\bf U}\right) &=& {\bf F}\left({\bf U}\right), \end{eqnarray}
where \({\bf F}\) is non-stiff and integrated in time explicitly, and \({\bf G}\) is stiff and integrated in time implicitly, and \({\bf U}\) represents the entire solution vector (state vector).
Note that \({\bf F}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time explicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).
This function computes \({\bf F}\left({\bf U}\right)\), given \({\bf U}\).
Note:
ts | Time integration object |
t | Current simulation time |
Y | State vector (input) |
F | The computed right-hand-side vector |
ctxt | Object of type PETScContext |
Definition at line 49 of file PetscRHSFunctionIMEX.cpp.
PetscErrorCode PetscIFunctionIMEX | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | Ydot, | ||
Vec | F, | ||
void * | ctxt | ||
) |
Compute the implicitly-treated part of the right-hand-side for the implicit-explicit (IMEX) time integration of the governing equations: The ODE, obtained after discretizing the governing PDE in space, is expressed as follows (for the purpose of IMEX time integration):
\begin{eqnarray} \frac {d{\bf U}}{dt} &=& {\bf F}\left({\bf U}\right) + {\bf G}\left({\bf U}\right), \\ \Rightarrow \dot{\bf U} - {\bf G}\left({\bf U}\right) &=& {\bf F}\left({\bf U}\right), \end{eqnarray}
where \({\bf F}\) is non-stiff and integrated in time explicitly, and \({\bf G}\) is stiff and integrated in time implicitly, and \({\bf U}\) represents the entire solution vector (state vector).
Note that \({\bf G}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time implicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).
This function computes the left-hand-side of the above equation:
\begin{equation} \mathcal{G}\left(\dot{\bf U},{\bf U},t\right) = \dot{\bf U} - {\bf G}\left({\bf U}\right) \end{equation}
given \(\dot{\bf U}\) and \({\bf U}\).
Notes:
ts | The time integration object |
t | Current solution time |
Y | State vector (input) |
Ydot | Time derivative of the state vector (input) |
F | The computed function vector |
ctxt | Object of type PETScContext |
Definition at line 53 of file PetscIFunctionIMEX.cpp.
PetscErrorCode PetscIFunctionImpl | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | Ydot, | ||
Vec | F, | ||
void * | ctxt | ||
) |
Compute the left-hand-side for the implicit time integration of the governing equations: The spatially discretized ODE can be expressed as
\begin{equation} \frac {d{\bf U}} {dt} = {\bf F}\left({\bf U}\right). \end{equation}
This function computes \(\dot{\bf U} - {\bf F}\left({\bf U}\right)\), given \({\bf U},\dot{\bf U}\).
Note:
ts | Time integration object |
t | Current simulation time |
Y | State vector (input) |
Ydot | Time derivative of the state vector (input) |
F | The computed right-hand-side vector |
ctxt | Object of type PETScContext |
Definition at line 37 of file PetscIFunctionImpl.cpp.
PetscErrorCode PetscIJacobianIMEX | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | Ydot, | ||
PetscReal | a, | ||
Mat | A, | ||
Mat | B, | ||
void * | ctxt | ||
) |
Compute the Jacobian for implicit-explicit (IMEX) 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) + {\bf G}\left({\bf U}\right) \Rightarrow \dot{\bf U} - {\bf G}\left({\bf U}\right) = {\bf F}\left({\bf U}\right), \end{equation}
where \({\bf F}\) is the spatially discretized non-stiff right-hand-side (treated explicitly), \({\bf G}\) is the spatially discretized stiff right-hand-side (treated implicitly), and \({\bf U}\) represents the entire solution vector. The Jacobian of the implicit part is thus given by:
\begin{equation} {\bf J} = \left[\alpha{\bf I} - \frac {\partial {\bf G}} {\partial {\bf U}} \right] \end{equation}
where \(\alpha\) is the shift coefficient (PETScContext::shift) of the time integration method.
Note that \({\bf G}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time implicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).
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 (PetscJacobianFunctionIMEX_JFNK() for nonlinear systems, and PetscJacobianFunctionIMEX_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 61 of file PetscIJacobianIMEX.cpp.
PetscErrorCode PetscJacobianFunctionIMEX_JFNK | ( | Mat | Jacobian, |
Vec | Y, | ||
Vec | F | ||
) |
Computes the action of the Jacobian on a vector: See documentation for PetscIJacobianIMEX() 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-explicit (IMEX) time integration of the ODE given by
\begin{equation} \frac {d {\bf U}} {dt} = {\bf F}\left({\bf U}\right) + {\bf G}\left({\bf U}\right) \Rightarrow \frac {d {\bf U}} {dt} - {\bf G}\left({\bf U}\right) = {\bf F}\left({\bf U}\right), \end{equation}
where \({\bf F}\) is the spatially discretized non-stiff right-hand-side (treated explicitly), \({\bf G}\) is the spatially discretized stiff right-hand-side (treated implicitly), we have that
\begin{equation} \mathcal{F}\left({\bf U}\right) \equiv \dot{\bf U} - {\bf G}\left({\bf U}\right) \Rightarrow {\bf J}\left({\bf U}\right) = \left[\alpha{\bf I} - \frac {\partial {\bf G}\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 G}\left({\bf U}_0+\epsilon{\bf y}\right)-{\bf G}\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.
Note that \({\bf G}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time implicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).
Notes:
Jacobian | Jacobian matrix |
Y | Input vector |
F | Output vector (Jacobian times input vector |
Definition at line 128 of file PetscIJacobianIMEX.cpp.
PetscErrorCode PetscJacobianFunctionIMEX_Linear | ( | Mat | Jacobian, |
Vec | Y, | ||
Vec | F | ||
) |
Computes the action of the Jacobian on a vector: See documentation for PetscIJacobianIMEX() 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-explicit (IMEX) time integration of the ODE given by
\begin{equation} \frac {d {\bf U}} {dt} = {\bf F}\left({\bf U}\right) + {\bf G}\left({\bf U}\right) \Rightarrow \frac {d {\bf U}} {dt} - {\bf G}\left({\bf U}\right) = {\bf F}\left({\bf U}\right), \end{equation}
where \({\bf F}\) is the spatially discretized non-stiff right-hand-side (treated explicitly), \({\bf G}\) is the spatially discretized stiff right-hand-side (treated implicitly) and is linear, we have that
\begin{equation} \mathcal{F}\left({\bf U}\right) \equiv \dot{\bf U} - {\bf G}\left({\bf U}\right) \Rightarrow {\bf J}\left({\bf U}\right) = \left[\alpha{\bf I} - \frac {\partial {\bf G}\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 G}\left({\bf U}_0+{\bf y}\right)-{\bf G}\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.
Note that \({\bf G}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time implicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).
Notes:
Jacobian | Jacobian matrix |
Y | Input vector |
F | Output vector (Jacobian times input vector |
Definition at line 290 of file PetscIJacobianIMEX.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.
int PetscGlobalDOF | ( | void * | c | ) |
Compute the global DOF index for all the grid points: The "global DOF index" is the component number (or block component number for HyPar::nvars > 1) of a grid point in the global solution vector. It is also the row number (or block row number) of the grid point in the global matrix representing, for example, the Jacobian of the right-hand-side.
PETScContext::globalDOF is an integer array with the same layout as the solution array HyPar::u (but with one component) containing the global DOF index for the corresponding grid points. It has the same number of ghost points as HyPar::u.
c | Object of type PETScContext |
Definition at line 69 of file PetscGlobalDOF.cpp.
int PetscCleanup | ( | void * | obj | ) |
Clean up allocations in the PETSc interface
obj | Object of type PETScContext |
Definition at line 12 of file PetscCleanup.cpp.
int PetscCreatePointList | ( | void * | obj | ) |
Create a list of computational points for each simulation domain: This is a list of all the grid points on which the PDE is solved. Thus, it is the total number of grid points minus the ghost points and blanked out points.
Note: this list is local, not global.
obj | Object of type PETScContext |
Definition at line 25 of file PetscCreatePointList.cpp.
PetscErrorCode PetscPreStage | ( | TS | ts, |
PetscReal | waqt | ||
) |
Function called before a stage in multi-stage time-integration methods
ts | Time integration object |
waqt | Current simulation time |
Definition at line 14 of file PetscPreStage.cpp.
PetscErrorCode PetscPostStage | ( | TS | ts, |
PetscReal | stagetime, | ||
PetscInt | stageindex, | ||
Vec * | Y | ||
) |
Function called after every stage in a multi-stage time-integration method
ts | Time integrator of PETSc type TS |
stagetime | Current stage time |
stageindex | Stage |
Y | Stage solutions (all stages) - be careful what you access |
Definition at line 16 of file PetscPostStage.cpp.
PetscErrorCode PetscPreTimeStep | ( | TS | ts | ) |
Function called before a time step
ts | Time integration object |
Definition at line 24 of file PetscPreTimeStep.cpp.
PetscErrorCode PetscPostTimeStep | ( | TS | ts | ) |
Function called after every time step
ts | Time integrator object |
Definition at line 21 of file PetscPostTimeStep.cpp.
int PetscComputePreconMatIMEX | ( | Mat | Pmat, |
Vec | Y, | ||
void * | ctxt | ||
) |
Compute and assemble the preconditioning matrix for the implicit-explicit (IMEX) time integration of the governing equations: Right now, it just calls PetscComputePreconMatImpl()
Pmat | Preconditioning matrix to construct |
Y | Solution vector |
ctxt | Application context |
Definition at line 17 of file PetscComputePreconMatIMEX.cpp.
int PetscComputePreconMatImpl | ( | Mat | Pmat, |
Vec | Y, | ||
void * | ctxt | ||
) |
Compute and assemble the preconditioning matrix for the 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 L}\left({\bf U}\right) \Rightarrow \frac {d{\bf U}}{dt} - {\bf L}\left({\bf U}\right) = 0, \end{equation}
where \({\bf L}\) 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 L}} {\partial {\bf U}} \right] \end{equation}
where \(\alpha\) is the shift coefficient (PETScContext::shift) of the time integration method.
Let \(\mathcal{D}_{\rm hyp}\) and \(\mathcal{D}_{\rm par}\) represent the spatial discretization methods for the hyperbolic and parabolic terms. Thus, the Jacobian can be written as follows:
\begin{equation} {\bf J} = \left[\alpha{\bf I} - \left(\mathcal{D}_{\rm hyp}\left\{\frac {\partial {\bf F}_{\rm hyp}} {\partial {\bf u}}\right\} + \mathcal{D}_{\rm par}\left\{\frac {\partial {\bf F}_{\rm par}} {\partial {\bf u}}\right\}\right)\right] \end{equation}
The preconditioning matrix is usually a close approximation of the actual Jacobian matrix, where the actual Jacobian may be too expensive to evaluate and assemble. In this function, the preconditioning matrix is the following approximation of the actual Jacobian:
\begin{equation} {\bf J}_p = \left[\alpha{\bf I} - \left(\mathcal{D}^{\left(l\right)}_{\rm hyp}\left\{\frac {\partial {\bf F}_{\rm hyp}} {\partial {\bf u}}\right\} + \mathcal{D}^{\left(l\right)}_{\rm par}\left\{\frac {\partial {\bf F}_{\rm par}} {\partial {\bf u}}\right\}\right) \right] \approx {\bf J}, \end{equation}
where \(\mathcal{D}^{\left(l\right)}_{\rm hyp,par}\) represent lower order discretizations of the hyperbolic and parabolic terms. The matrix \({\bf J}_p\) is provided to the preconditioner.
Note that the specific physics model provides the following functions:
Currently, this function doesn't include the source term.
Pmat | Preconditioning matrix to construct |
Y | Solution vector |
ctxt | Application context |
Definition at line 59 of file PetscComputePreconMatImpl.cpp.
int PetscJacobianMatNonzeroEntriesImpl | ( | Mat | Amat, |
int | width, | ||
void * | ctxt | ||
) |
Set non-zero entries of the Jacobian matrix for the 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 L}\left({\bf U}\right) \Rightarrow \frac {d{\bf U}}{dt} - {\bf L}\left({\bf U}\right) = 0, \end{equation}
where \({\bf L}\) 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 L}} {\partial {\bf U}} \right] \end{equation}
where \(\alpha\) is the shift coefficient (PETScContext::shift) of the time integration method.
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.
Amat | Matrix |
width | Stencil width |
ctxt | Application context |
Definition at line 38 of file PetscJacobianMatNonzeroEntriesImpl.cpp.
PetscErrorCode PetscTimeError | ( | TS | ts | ) |
Function to compute any error estimates, if available
Compute the norms of the error estimate, if the PETSc time integrator has it (for example TSGLEE)
ts | Time integrator object of PETSc type TS |
Definition at line 20 of file PetscError.cpp.
PetscErrorCode PetscSetInitialGuessROM | ( | SNES | snes, |
Vec | X, | ||
void * | ctxt | ||
) |
Function to compute initial guess of nonlinear solve using libROM
Compute the initial guess for a nonlinear solve using a trained libROM reduced-order model.
snes | Nonlinear solver object (see PETSc SNES) |
X | Initial guess vector |
ctxt | Application context |
Definition at line 18 of file PetscSetInitialGuessROM.cpp.