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 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);

 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: }

 50: /* The following are different routines used to modify the preconditioner */

 52: /*@

 54:   KSPFGMRESModifyPCNoChange - this is the default used by fgmres - it doesn't change the preconditioner.

 56:   Input Parameters:
 57: +    ksp - the ksp context being used.
 58: .    total_its     - the total number of FGMRES iterations that have occurred.
 59: .    loc_its       - the number of FGMRES iterations since last restart.
 60:                     a restart (so number of Krylov directions to be computed)
 61: .    res_norm      - the current residual norm.
 62: -    dummy         - context variable, unused in this routine

 64:    Level: intermediate

 66:    Contributed by Allison Baker

 68: You can use this as a template!

 70: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()

 72: @*/
 73: PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
 74: {
 76:   return(0);
 77: }

 79: /*@

 81:  KSPFGMRESModifyPCKSP - modifies the attributes of the
 82:      GMRES preconditioner.  It serves as an example (not as something
 83:      useful!)

 85:   Input Parameters:
 86: +    ksp - the ksp context being used.
 87: .    total_its     - the total number of FGMRES iterations that have occurred.
 88: .    loc_its       - the number of FGMRES iterations since last restart.
 89: .    res_norm      - the current residual norm.
 90: -    dummy         - context, not used here

 92:    Level: intermediate

 94:    Contributed by Allison Baker

 96:  This could be used as a template!

 98: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()

100: @*/
101: PetscErrorCode  KSPFGMRESModifyPCKSP(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
102: {
103:   PC             pc;
105:   PetscInt       maxits;
106:   KSP            sub_ksp;
107:   PetscReal      rtol,abstol,dtol;
108:   PetscBool      isksp;

111:   KSPGetPC(ksp,&pc);

113:   PetscObjectTypeCompare((PetscObject)pc,PCKSP,&isksp);
114:   if (isksp) {
115:     PCKSPGetKSP(pc,&sub_ksp);

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

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

122:     KSPGetTolerances(sub_ksp,&rtol,&abstol,&dtol,&maxits);
123:     if (!loc_its) rtol = .1;
124:     else rtol *= .9;
125:     KSPSetTolerances(sub_ksp,rtol,abstol,dtol,maxits);
126:   }
127:   return(0);
128: }