Actual source code: modpcf.c


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

  4: /*@C
  5:    KSPFGMRESSetModifyPC - Sets the routine used by `KSPFGMRES` to modify the preconditioner. [](sec_flexibleksp)

  7:    Logically Collective

  9:    Input Parameters:
 10: +  ksp - iterative context obtained from `KSPCreate()`
 11: .  fcn - modifypc function
 12: .  ctx - optional context
 13: -  d - optional context destroy routine

 15:    Calling Sequence of function:
 16:     PetscErrorCode fcn(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void*ctx);
 17: +    ksp - the ksp context being used.
 18: .    total_its     - the total number of FGMRES iterations that have occurred.
 19: .    loc_its       - the number of FGMRES iterations since last restart.
 20: .    res_norm      - the current residual norm.
 21: -    ctx           - optional context variable

 23:    Options Database Keys:
 24: +   -ksp_fgmres_modifypcnochange - do not change the `PC`
 25: -   -ksp_fgmres_modifypcksp - changes the inner KSP solver tolerances

 27:    Level: intermediate

 29:    Note:
 30:    Several modifypc routines are predefined, including  `KSPFGMRESModifyPCNoChange()`, and  `KSPFGMRESModifyPCKSP()`

 32: .seealso: [](chapter_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFGMRESModifyPCNoChange()`, `KSPFGMRESModifyPCKSP()`
 33: @*/
 34: PetscErrorCode KSPFGMRESSetModifyPC(KSP ksp, PetscErrorCode (*fcn)(KSP, PetscInt, PetscInt, PetscReal, void *), void *ctx, PetscErrorCode (*d)(void *))
 35: {
 36:   PetscFunctionBegin;
 38:   PetscTryMethod(ksp, "KSPFGMRESSetModifyPC_C", (KSP, PetscErrorCode(*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)), (ksp, fcn, ctx, d));
 39:   PetscFunctionReturn(PETSC_SUCCESS);
 40: }

 42: /*@
 43:   KSPFGMRESModifyPCNoChange - this is the default used by `KSPFMGMRES` - it doesn't change the preconditioner. [](sec_flexibleksp)

 45:   Input Parameters:
 46: +    ksp - the ksp context being used.
 47: .    total_its     - the total number of `KSPFGMRES` iterations that have occurred.
 48: .    loc_its       - the number of `KSPFGMRES` iterations since last restart.
 49:                     a restart (so number of Krylov directions to be computed)
 50: .    res_norm      - the current residual norm.
 51: -    dummy         - context variable, unused in this routine

 53:    Level: intermediate

 55: .seealso: [](chapter_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFGMRESSetModifyPC()`, `KSPFGMRESModifyPCKSP()`
 56: @*/
 57: PetscErrorCode KSPFGMRESModifyPCNoChange(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, void *dummy)
 58: {
 59:   PetscFunctionBegin;
 60:   PetscFunctionReturn(PETSC_SUCCESS);
 61: }

 63: /*@
 64:      KSPFGMRESModifyPCKSP - modifies the attributes of the `KSPFGMRES` preconditioner. [](sec_flexibleksp). It serves as an example (not as something useful in practice)

 66:   Input Parameters:
 67: +    ksp - the ksp context being used.
 68: .    total_its     - the total number of `KSPFGMRES` iterations that have occurred.
 69: .    loc_its       - the number of `KSPFGMRES` iterations since last restart.
 70: .    res_norm      - the current residual norm.
 71: -    dummy         - context, not used here

 73:    Level: intermediate

 75:    Note:
 76:     You can use this as a template for writing a custom monification callback

 78: .seealso: [](chapter_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFGMRESSetModifyPC()`, `KSPFGMRESModifyPCKSP()`
 79: @*/
 80: PetscErrorCode KSPFGMRESModifyPCKSP(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, void *dummy)
 81: {
 82:   PC        pc;
 83:   PetscInt  maxits;
 84:   KSP       sub_ksp;
 85:   PetscReal rtol, abstol, dtol;
 86:   PetscBool isksp;

 88:   PetscFunctionBegin;
 89:   PetscCall(KSPGetPC(ksp, &pc));

 91:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCKSP, &isksp));
 92:   if (isksp) {
 93:     PetscCall(PCKSPGetKSP(pc, &sub_ksp));

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

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

100:     PetscCall(KSPGetTolerances(sub_ksp, &rtol, &abstol, &dtol, &maxits));
101:     if (!loc_its) rtol = .1;
102:     else rtol *= .9;
103:     PetscCall(KSPSetTolerances(sub_ksp, rtol, abstol, dtol, maxits));
104:   }
105:   PetscFunctionReturn(PETSC_SUCCESS);
106: }