Actual source code: pcsles.c

  1: /*$Id: pcsles.c,v 1.39 2001/04/10 19:36:17 bsmith Exp $*/


 4:  #include src/sles/pc/pcimpl.h
 5:  #include petscsles.h

  7: typedef struct {
  8:   PetscTruth use_true_matrix;       /* use mat rather than pmat in inner linear solve */
  9:   SLES       sles;
 10:   int        its;                   /* total number of iterations SLES uses */
 11: } PC_SLES;

 13: #undef __FUNCT__  
 15: static int PCApply_SLES(PC pc,Vec x,Vec y)
 16: {
 17:   int     ierr,its;
 18:   PC_SLES *jac = (PC_SLES*)pc->data;

 21:   ierr      = SLESSolve(jac->sles,x,y,&its);
 22:   jac->its += its;
 23:   return(0);
 24: }

 26: #undef __FUNCT__  
 28: static int PCApplyTranspose_SLES(PC pc,Vec x,Vec y)
 29: {
 30:   int     ierr,its;
 31:   PC_SLES *jac = (PC_SLES*)pc->data;

 34:   ierr      = SLESSolveTranspose(jac->sles,x,y,&its);
 35:   jac->its += its;
 36:   return(0);
 37: }

 39: #undef __FUNCT__  
 41: static int PCSetUp_SLES(PC pc)
 42: {
 43:   int     ierr;
 44:   PC_SLES *jac = (PC_SLES*)pc->data;
 45:   Mat     mat;

 48:   SLESSetFromOptions(jac->sles);
 49:   if (jac->use_true_matrix) mat = pc->mat;
 50:   else                      mat = pc->pmat;

 52:   SLESSetOperators(jac->sles,mat,pc->pmat,pc->flag);
 53:   SLESSetUp(jac->sles,pc->vec,pc->vec);
 54:   return(0);
 55: }

 57: /* Default destroy, if it has never been setup */
 58: #undef __FUNCT__  
 60: static int PCDestroy_SLES(PC pc)
 61: {
 62:   PC_SLES *jac = (PC_SLES*)pc->data;
 63:   int     ierr;

 66:   SLESDestroy(jac->sles);
 67:   PetscFree(jac);
 68:   return(0);
 69: }

 71: #undef __FUNCT__  
 73: static int PCView_SLES(PC pc,PetscViewer viewer)
 74: {
 75:   PC_SLES    *jac = (PC_SLES*)pc->data;
 76:   int        ierr;
 77:   PetscTruth isascii;

 80:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
 81:   if (isascii) {
 82:     if (jac->use_true_matrix) {
 83:       PetscViewerASCIIPrintf(viewer,"Using true matrix (not preconditioner matrix) on inner solven");
 84:     }
 85:     PetscViewerASCIIPrintf(viewer,"KSP and PC on SLES preconditioner follown");
 86:     PetscViewerASCIIPrintf(viewer,"---------------------------------n");
 87:   } else {
 88:     SETERRQ1(1,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
 89:   }
 90:   PetscViewerASCIIPushTab(viewer);
 91:   SLESView(jac->sles,viewer);
 92:   PetscViewerASCIIPopTab(viewer);
 93:   if (isascii) {
 94:     PetscViewerASCIIPrintf(viewer,"---------------------------------n");
 95:   }
 96:   return(0);
 97: }

 99: #undef __FUNCT__  
101: static int PCSetFromOptions_SLES(PC pc){
102:   int        ierr;
103:   PetscTruth flg;

106:   PetscOptionsHead("SLES preconditioner options");
107:     PetscOptionsName("-pc_sles_true","Use true matrix to define inner linear system, not preconditioner matrix","PCSLESSetUseTrue",&flg);
108:     if (flg) {
109:       PCSLESSetUseTrue(pc);
110:     }
111:   PetscOptionsTail();
112:   return(0);
113: }

115: /* ----------------------------------------------------------------------------------*/

117: EXTERN_C_BEGIN
118: #undef __FUNCT__  
120: int PCSLESSetUseTrue_SLES(PC pc)
121: {
122:   PC_SLES   *jac;

125:   jac                  = (PC_SLES*)pc->data;
126:   jac->use_true_matrix = PETSC_TRUE;
127:   return(0);
128: }
129: EXTERN_C_END

131: EXTERN_C_BEGIN
132: #undef __FUNCT__  
134: int PCSLESGetSLES_SLES(PC pc,SLES *sles)
135: {
136:   PC_SLES   *jac;

139:   jac          = (PC_SLES*)pc->data;
140:   *sles        = jac->sles;
141:   return(0);
142: }
143: EXTERN_C_END

145: #undef __FUNCT__  
147: /*@
148:    PCSLESSetUseTrue - Sets a flag to indicate that the true matrix (rather than
149:    the matrix used to define the preconditioner) is used to compute the
150:    residual inside the inner solve.

152:    Collective on PC

154:    Input Parameters:
155: .  pc - the preconditioner context

157:    Options Database Key:
158: .  -pc_sles_true - Activates PCSLESSetUseTrue()

160:    Note:
161:    For the common case in which the preconditioning and linear 
162:    system matrices are identical, this routine is unnecessary.

164:    Level: advanced

166: .keywords:  PC, SLES, set, true, local, flag

168: .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal()
169: @*/
170: int PCSLESSetUseTrue(PC pc)
171: {
172:   int ierr,(*f)(PC);

176:   PetscObjectQueryFunction((PetscObject)pc,"PCSLESSetUseTrue_C",(void (**)(void))&f);
177:   if (f) {
178:     (*f)(pc);
179:   }
180:   return(0);
181: }

183: #undef __FUNCT__  
185: /*@C
186:    PCSLESGetSLES - Gets the SLES context for a SLES PC.

188:    Not Collective but SLES returned is parallel if PC was parallel

190:    Input Parameter:
191: .  pc - the preconditioner context

193:    Output Parameters:
194: .  sles - the PC solver

196:    Notes:
197:    You must call SLESSetUp() before calling PCSLESGetSLES().

199:    Level: advanced

201: .keywords:  PC, SLES, get, context
202: @*/
203: int PCSLESGetSLES(PC pc,SLES *sles)
204: {
205:   int ierr,(*f)(PC,SLES*);

209:   if (!pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SLESSetUp first");
210:   PetscObjectQueryFunction((PetscObject)pc,"PCSLESGetSLES_C",(void (**)(void))&f);
211:   if (f) {
212:     (*f)(pc,sles);
213:   }
214:   return(0);
215: }

217: /* ----------------------------------------------------------------------------------*/

219: /*MC
220:      PCSLES -    Defines a preconditioner that can consist of any SLES solver.
221:                  This allows, for example, embedding a Krylov method inside a preconditioner.

223:    Options Database Key:
224: .     -pc_sles_true - use the matrix that defines the linear system as the matrix for the
225:                     inner solver, otherwise by default it uses the matrix used to construct
226:                     the preconditioner (see PCSetOperators())

228:    Level: intermediate

230:    Concepts: inner iteration

232:    Notes: Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
233:           the incorrect answer) unless you use KSPFGMRES as the other Krylov method


236: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
237:            PCSHELL, PCCOMPOSITE, PCSLESUseTrue(), PCSLESGetSLES()

239: M*/

241: EXTERN_C_BEGIN
242: #undef __FUNCT__  
244: int PCCreate_SLES(PC pc)
245: {
246:   int       ierr;
247:   char      *prefix;
248:   PC_SLES   *jac;

251:   PetscNew(PC_SLES,&jac);
252:   PetscLogObjectMemory(pc,sizeof(PC_SLES));
253:   pc->ops->apply              = PCApply_SLES;
254:   pc->ops->applytranspose     = PCApplyTranspose_SLES;
255:   pc->ops->setup              = PCSetUp_SLES;
256:   pc->ops->destroy            = PCDestroy_SLES;
257:   pc->ops->setfromoptions     = PCSetFromOptions_SLES;
258:   pc->ops->view               = PCView_SLES;
259:   pc->ops->applyrichardson    = 0;

261:   pc->data               = (void*)jac;
262:   ierr                   = SLESCreate(pc->comm,&jac->sles);

264:   PCGetOptionsPrefix(pc,&prefix);
265:   SLESSetOptionsPrefix(jac->sles,prefix);
266:   SLESAppendOptionsPrefix(jac->sles,"sles_");
267:   jac->use_true_matrix = PETSC_FALSE;
268:   jac->its             = 0;

270:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSLESSetUseTrue_C","PCSLESSetUseTrue_SLES",
271:                     PCSLESSetUseTrue_SLES);
272:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSLESGetSLES_C","PCSLESGetSLES_SLES",
273:                     PCSLESGetSLES_SLES);

275:   return(0);
276: }
277: EXTERN_C_END