Actual source code: modpcf.c

  1: /* $Id: modpcf.c,v 1.14 2001/04/10 19:36:35 bsmith Exp $*/

 3:  #include petscsles.h
  4: #undef __FUNCT__  
  6: /*@C
  7:    KSPFGMRESSetModifyPC - Sets the routine used by FGMRES to modify the preconditioner.

  9:    Collective on KSP

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

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

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

 26:    Options Database Keys:
 27:    -ksp_fgmres_modifypcnochange
 28:    -ksp_fgmres_modifypcsles

 30:    Level: intermediate

 32:    Contributed by Allison Baker

 34:    Notes:
 35:    Several modifypc routines are predefined, including
 36:     KSPFGMRESModifyPCNoChange()
 37:     KSPFGMRESModifyPCSLES()

 39: .seealso: KSPFGMRESModifyPCNoChange(), KSPFGMRESModifyPCSLES()

 41: @*/
 42: int KSPFGMRESSetModifyPC(KSP ksp,int (*fcn)(KSP,int,int,PetscReal,void*),void* ctx,int (*d)(void*))
 43: {
 44:   int ierr,(*f)(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int (*)(void*));

 48:   PetscObjectQueryFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",(void (**)(void))&f);
 49:   if (f) {
 50:     (*f)(ksp,fcn,ctx,d);
 51:   }
 52:   return(0);
 53: }



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

 59: #undef __FUNCT__  
 61: /*@C

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

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

 73:    Level: intermediate

 75:    Contributed by Allison Baker

 77: You can use this as a template!

 79: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCSLES()

 81: @*/
 82: int KSPFGMRESModifyPCNoChange(KSP ksp,int total_its,int loc_its,PetscReal res_norm,void* dummy)
 83: {

 86:   return(0);
 87: }

 89: #undef __FUNCT__  
 91: /*@C

 93:  KSPFGMRESModifyPCSLES - modifies the attributes of the
 94:      GMRES preconditioner.  It serves as an example (not as something 
 95:      useful!) 

 97:   Input Parameters:
 98: +    ksp - the ksp context being used.
 99: .    total_its     - the total number of FGMRES iterations that have occurred.    
100: .    loc_its       - the number of FGMRES iterations since last restart.
101: .    res_norm      - the current residual norm.
102: -    dummy         - context, not used here

104:    Level: intermediate

106:    Contributed by Allison Baker

108:  This could be used as a template!

110: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCSLES()

112: @*/
113: int KSPFGMRESModifyPCSLES(KSP ksp,int total_its,int loc_its,PetscReal res_norm,void *dummy)
114: {
115:   PC         pc;
116:   int        ierr,maxits;
117:   SLES       sub_sles;
118:   KSP        sub_ksp;
119:   PetscReal  rtol,atol,dtol;
120:   PetscTruth issles;


124:   KSPGetPC(ksp,&pc);

126:   PetscTypeCompare((PetscObject)pc,PCSLES,&issles);
127:   if (issles) {
128:     PCSLESGetSLES(pc,&sub_sles);
129:     SLESGetKSP(sub_sles,&sub_ksp);
130: 
131:     /* note that at this point you could check the type of KSP with KSPGetType() */

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

136:     KSPGetTolerances(sub_ksp,&rtol,&atol,&dtol,&maxits);
137:     if (loc_its == 0) {
138:       rtol = .1;
139:     } else {
140:       rtol *= .9;
141:     }
142:     KSPSetTolerances(sub_ksp,rtol,atol,dtol,maxits);
143:   }
144:   return(0);
145: }