Actual source code: modpcf.c

  1: #include <petsc/private/kspimpl.h>

  3: /*@C
  4:   KSPFlexibleSetModifyPC - Sets the routine used by flexible `KSP` methods to modify the preconditioner. [](sec_flexibleksp)

  6:   Logically Collective

  8:   Input Parameters:
  9: + ksp     - iterative context obtained from `KSPCreate()`
 10: . fcn     - function to modify the `PC`, see `KSPFlexibleModifyPCFn`
 11: . ctx     - optional context
 12: - destroy - optional context destroy routine

 14:   Level: intermediate

 16:   Note:
 17:   Several `fcn` routines are predefined, including `KSPFlexibleModifyPCNoChange()` and `KSPFlexibleModifyPCKSP()`

 19: .seealso: [](ch_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFCG`, `KSPPIPEFCG`, `KSPGCR`, `KSPPIPEGCR`, `KSPFlexibleModifyPCFn`, `KSPFlexibleModifyPCNoChange()`, `KSPFlexibleModifyPCKSP()`
 20: @*/
 21: PetscErrorCode KSPFlexibleSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fcn, PetscCtx ctx, PetscCtxDestroyFn *destroy)
 22: {
 23:   PetscFunctionBegin;
 25:   PetscTryMethod(ksp, "KSPFlexibleSetModifyPC_C", (KSP, KSPFlexibleModifyPCFn *, PetscCtx, PetscCtxDestroyFn *), (ksp, fcn, ctx, destroy));
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: /*@C
 30:   KSPFlexibleModifyPCNoChange - this is the default used by the flexible Krylov methods - it doesn't change the preconditioner. [](sec_flexibleksp)

 32:   Input Parameters:
 33: + ksp       - the ksp context being used.
 34: . total_its - the total number of `KSP` iterations that have occurred.
 35: . loc_its   - the number of `KSP` iterations since last restart.
 36: . res_norm  - the current residual norm.
 37: - ctx       - context variable, unused in this routine

 39:   Level: intermediate

 41:   Note:
 42:   See `KSPFlexibleModifyPCKSP()` for a template for providing your own modification routines.

 44: .seealso: [](ch_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFCG`, `KSPPIPEFCG`, `KSPGCR`, `KSPPIPEGCR`, `KSPFlexibleModifyPCFn`, `KSPFlexibleSetModifyPC()`, `KSPFlexibleModifyPCKSP()`
 45: @*/
 46: PetscErrorCode KSPFlexibleModifyPCNoChange(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx)
 47: {
 48:   PetscFunctionBegin;
 49:   PetscFunctionReturn(PETSC_SUCCESS);
 50: }

 52: /*@C
 53:   KSPFlexibleModifyPCKSP - modifies the attributes of the `PCKSP` preconditioner, see [](sec_flexibleksp).

 55:   Input Parameters:
 56: + ksp       - the ksp context being used.
 57: . total_its - the total number of `KSP` iterations that have occurred.
 58: . loc_its   - the number of `KSP` iterations since last restart.
 59: . res_norm  - the current residual norm.
 60: - ctx       - context, unused in this routine

 62:   Level: intermediate

 64:   Notes:
 65:   You can use this as a template for writing a custom modification callback. See the source code for the change it makes.

 67:   You can provide this to a `KSP` with `KSPFlexibleSetModifyPC()`

 69: .seealso: [](ch_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFCG`, `KSPPIPEFCG`, `KSPGCR`, `KSPPIPEGCR`, `KSPFlexibleModifyPCFn`, `KSPFlexibleSetModifyPC()`
 70: @*/
 71: PetscErrorCode KSPFlexibleModifyPCKSP(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx)
 72: {
 73:   PC        pc;
 74:   PetscInt  maxits;
 75:   KSP       sub_ksp;
 76:   PetscReal rtol, abstol, dtol;
 77:   PetscBool isksp;

 79:   PetscFunctionBegin;
 80:   PetscCall(KSPGetPC(ksp, &pc));

 82:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCKSP, &isksp));
 83:   if (isksp) {
 84:     PetscCall(PCKSPGetKSP(pc, &sub_ksp));

 86:     /* note that at this point you could check the type of KSP with KSPGetType() */

 88:     /* Now we can use functions such as KSPGMRESSetRestart() or
 89:       KSPGMRESSetOrthogonalization() or KSPSetTolerances() */

 91:     PetscCall(KSPGetTolerances(sub_ksp, &rtol, &abstol, &dtol, &maxits));
 92:     if (!loc_its) rtol = .1;
 93:     else rtol *= .9;
 94:     PetscCall(KSPSetTolerances(sub_ksp, rtol, abstol, dtol, maxits));
 95:   }
 96:   PetscFunctionReturn(PETSC_SUCCESS);
 97: }