petsc-3.14.6 2021-03-30
PCGalerkinSetComputeSubmatrix
Provide a routine that will be called to compute the Galerkin submatrix
Synopsis
#include "petscksp.h"
PetscErrorCode PCGalerkinSetComputeSubmatrix(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx)
Logically Collective
Input Parameter
| pc | - the preconditioner context
|
| computeAsub | - routine that computes the submatrix from the global matrix
|
| ctx | - context used by the routine, or NULL
|
Calling sequence of computeAsub
computeAsub(PC pc,Mat A, Mat Ap, Mat *cAP,void *ctx);
| PC | - the Galerkin PC
|
| A | - the matrix in the Galerkin PC
|
| Ap | - the computed submatrix from any previous computation, if NULL it has not previously been computed
|
| cAp | - the submatrix computed by this routine
|
| ctx | - optional user-defined function context
|
Notes
Instead of providing this routine you can call PCGalerkinGetKSP() and then KSPSetOperators() to provide the submatrix,
but that will not work for multiple KSPSolves with different matrices unless you call it for each solve.
This routine is called each time the outer matrix is changed. In the first call the Ap argument is NULL and the routine should create the
matrix and computes its values in cAp. On each subsequent call the routine should up the Ap matrix.
Developer Notes
If the user does not call this routine nor call PCGalerkinGetKSP() and KSPSetOperators() then PCGalerkin could
could automatically compute the submatrix via calls to MatGalerkin() or MatRARt()
See Also
PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP()
Level
Intermediate
Location
src/ksp/pc/impls/galerkin/galerkin.c
Implementations
PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (*computeAsub) in src/ksp/pc/impls/galerkin/galerkin.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages