Actual source code: modpcf.c

petsc-3.3-p7 2013-05-11
  2: #include <petscksp.h>  /*I "petscksp.h" I*/
  5: /*@C
  6:    KSPFGMRESSetModifyPC - Sets the routine used by FGMRES to modify the preconditioner.

  8:    Logically Collective on KSP

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

 16:    Calling Sequence of function:
 17:     int fcn(KSP ksp,int total_its,int loc_its,PetscReal res_norm,void*ctx);

 19:     ksp - the ksp context being used.
 20:     total_its     - the total number of FGMRES iterations that have occurred.    
 21:     loc_its       - the number of FGMRES iterations since last restart.
 22:     res_norm      - the current residual norm.
 23:     ctx           - optional context variable

 25:    Options Database Keys:
 26:    -ksp_fgmres_modifypcnochange
 27:    -ksp_fgmres_modifypcksp

 29:    Level: intermediate

 31:    Contributed by Allison Baker

 33:    Notes:
 34:    Several modifypc routines are predefined, including
 35:     KSPFGMRESModifyPCNoChange()
 36:     KSPFGMRESModifyPCKSP()

 38: .seealso: KSPFGMRESModifyPCNoChange(), KSPFGMRESModifyPCKSP()

 40: @*/
 41: PetscErrorCode  KSPFGMRESSetModifyPC(KSP ksp,PetscErrorCode (*fcn)(KSP,PetscInt,PetscInt,PetscReal,void*),void* ctx,PetscErrorCode (*d)(void*))
 42: {

 47:   PetscTryMethod(ksp,"KSPFGMRESSetModifyPC_C",(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void*)),(ksp,fcn,ctx,d));
 48:   return(0);
 49: }


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

 56: /*@

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

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

 68:    Level: intermediate

 70:    Contributed by Allison Baker

 72: You can use this as a template!

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

 76: @*/
 77: PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void* dummy)
 78: {

 81:   return(0);
 82: }

 86: /*@

 88:  KSPFGMRESModifyPCKSP - modifies the attributes of the
 89:      GMRES preconditioner.  It serves as an example (not as something 
 90:      useful!) 

 92:   Input Parameters:
 93: +    ksp - the ksp context being used.
 94: .    total_its     - the total number of FGMRES iterations that have occurred.    
 95: .    loc_its       - the number of FGMRES iterations since last restart.
 96: .    res_norm      - the current residual norm.
 97: -    dummy         - context, not used here

 99:    Level: intermediate

101:    Contributed by Allison Baker

103:  This could be used as a template!

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

107: @*/
108: PetscErrorCode  KSPFGMRESModifyPCKSP(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
109: {
110:   PC             pc;
112:   PetscInt       maxits;
113:   KSP            sub_ksp;
114:   PetscReal      rtol,abstol,dtol;
115:   PetscBool      isksp;

118:   KSPGetPC(ksp,&pc);

120:   PetscObjectTypeCompare((PetscObject)pc,PCKSP,&isksp);
121:   if (isksp) {
122:     PCKSPGetKSP(pc,&sub_ksp);
123: 
124:     /* note that at this point you could check the type of KSP with KSPGetType() */

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

129:     KSPGetTolerances(sub_ksp,&rtol,&abstol,&dtol,&maxits);
130:     if (!loc_its) {
131:       rtol = .1;
132:     } else {
133:       rtol *= .9;
134:     }
135:     KSPSetTolerances(sub_ksp,rtol,abstol,dtol,maxits);
136:   }
137:   return(0);
138: }