Actual source code: modpcf.c
  2: #include <petsc/private/kspimpl.h>
  4: /*@C
  5:    KSPFGMRESSetModifyPC - Sets the routine used by FGMRES to modify the preconditioner.
  7:    Logically Collective on ksp
  9:    Input Parameters:
 10: +  ksp - iterative context obtained from KSPCreate
 11: .  fcn - modifypc function
 12: .  ctx - optional contex
 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);
 18:     ksp - the ksp context being used.
 19:     total_its     - the total number of FGMRES iterations that have occurred.
 20:     loc_its       - the number of FGMRES iterations since last restart.
 21:     res_norm      - the current residual norm.
 22:     ctx           - optional context variable
 24:    Options Database Keys:
 25:    -ksp_fgmres_modifypcnochange
 26:    -ksp_fgmres_modifypcksp
 28:    Level: intermediate
 30:    Contributed by Allison Baker
 32:    Notes:
 33:    Several modifypc routines are predefined, including
 34:     KSPFGMRESModifyPCNoChange()
 35:     KSPFGMRESModifyPCKSP()
 37: .seealso: KSPFGMRESModifyPCNoChange(), KSPFGMRESModifyPCKSP()
 39: @*/
 40: PetscErrorCode  KSPFGMRESSetModifyPC(KSP ksp,PetscErrorCode (*fcn)(KSP,PetscInt,PetscInt,PetscReal,void*),void *ctx,PetscErrorCode (*d)(void*))
 41: {
 46:   PetscTryMethod(ksp,"KSPFGMRESSetModifyPC_C",(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void*)),(ksp,fcn,ctx,d));
 47:   return(0);
 48: }
 51: /* The following are different routines used to modify the preconditioner */
 53: /*@
 55:   KSPFGMRESModifyPCNoChange - this is the default used by fgmres - it doesn't change the preconditioner.
 57:   Input Parameters:
 58: +    ksp - the ksp context being used.
 59: .    total_its     - the total number of FGMRES iterations that have occurred.
 60: .    loc_its       - the number of FGMRES iterations since last restart.
 61:                     a restart (so number of Krylov directions to be computed)
 62: .    res_norm      - the current residual norm.
 63: -    dummy         - context variable, unused in this routine
 65:    Level: intermediate
 67:    Contributed by Allison Baker
 69: You can use this as a template!
 71: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()
 73: @*/
 74: PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
 75: {
 77:   return(0);
 78: }
 80: /*@
 82:  KSPFGMRESModifyPCKSP - modifies the attributes of the
 83:      GMRES preconditioner.  It serves as an example (not as something
 84:      useful!)
 86:   Input Parameters:
 87: +    ksp - the ksp context being used.
 88: .    total_its     - the total number of FGMRES iterations that have occurred.
 89: .    loc_its       - the number of FGMRES iterations since last restart.
 90: .    res_norm      - the current residual norm.
 91: -    dummy         - context, not used here
 93:    Level: intermediate
 95:    Contributed by Allison Baker
 97:  This could be used as a template!
 99: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()
101: @*/
102: PetscErrorCode  KSPFGMRESModifyPCKSP(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
103: {
104:   PC             pc;
106:   PetscInt       maxits;
107:   KSP            sub_ksp;
108:   PetscReal      rtol,abstol,dtol;
109:   PetscBool      isksp;
112:   KSPGetPC(ksp,&pc);
114:   PetscObjectTypeCompare((PetscObject)pc,PCKSP,&isksp);
115:   if (isksp) {
116:     PCKSPGetKSP(pc,&sub_ksp);
118:     /* note that at this point you could check the type of KSP with KSPGetType() */
120:     /* Now we can use functions such as KSPGMRESSetRestart() or
121:       KSPGMRESSetOrthogonalization() or KSPSetTolerances() */
123:     KSPGetTolerances(sub_ksp,&rtol,&abstol,&dtol,&maxits);
124:     if (!loc_its) rtol = .1;
125:     else rtol *= .9;
126:     KSPSetTolerances(sub_ksp,rtol,abstol,dtol,maxits);
127:   }
128:   return(0);
129: }