Actual source code: fieldsplit.c

petsc-master 2016-07-27
Report Typos and Errors
  2: #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
  3: #include <petsc/private/kspimpl.h>    /*  This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
  4: #include <petscdm.h>

  6: const char *const PCFieldSplitSchurPreTypes[] = {"SELF","SELFP","A11","USER","FULL","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
  7: const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};

  9: PetscLogEvent KSP_Solve_FS_0,KSP_Solve_FS_1,KSP_Solve_FS_S,KSP_Solve_FS_U,KSP_Solve_FS_L,KSP_Solve_FS_2,KSP_Solve_FS_3,KSP_Solve_FS_4;

 11: typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
 12: struct _PC_FieldSplitLink {
 13:   KSP               ksp;
 14:   Vec               x,y,z;
 15:   char              *splitname;
 16:   PetscInt          nfields;
 17:   PetscInt          *fields,*fields_col;
 18:   VecScatter        sctx;
 19:   IS                is,is_col,is_orig;
 20:   PC_FieldSplitLink next,previous;
 21:   PetscLogEvent     event;
 22: };

 24: typedef struct {
 25:   PCCompositeType type;
 26:   PetscBool       defaultsplit;                    /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
 27:   PetscBool       splitdefined;                    /* Flag is set after the splits have been defined, to prevent more splits from being added */
 28:   PetscInt        bs;                              /* Block size for IS and Mat structures */
 29:   PetscInt        nsplits;                         /* Number of field divisions defined */
 30:   Vec             *x,*y,w1,w2;
 31:   Mat             *mat;                            /* The diagonal block for each split */
 32:   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
 33:   Mat             *Afield;                         /* The rows of the matrix associated with each split */
 34:   PetscBool       issetup;

 36:   /* Only used when Schur complement preconditioning is used */
 37:   Mat                       B;                     /* The (0,1) block */
 38:   Mat                       C;                     /* The (1,0) block */
 39:   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
 40:   Mat                       schurp;                /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
 41:   Mat                       schur_user;            /* User-provided preconditioning matrix for the Schur complement */
 42:   PCFieldSplitSchurPreType  schurpre;              /* Determines which preconditioning matrix is used for the Schur complement */
 43:   PCFieldSplitSchurFactType schurfactorization;
 44:   KSP                       kspschur;              /* The solver for S */
 45:   KSP                       kspupper;              /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
 46:   PC_FieldSplitLink         head;
 47:   PetscBool                 reset;                  /* indicates PCReset() has been last called on this object, hack */
 48:   PetscBool                 isrestrict;             /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
 49:   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
 50:   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
 51:   PetscBool                 diag_use_amat;          /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
 52:   PetscBool                 offdiag_use_amat;       /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
 53: } PC_FieldSplit;

 55: /*
 56:     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
 57:    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
 58:    PC you could change this.
 59: */

 61: /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
 62: * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
 63: static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
 64: {
 65:   switch (jac->schurpre) {
 66:   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
 67:   case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
 68:   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
 69:   case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
 70:   case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
 71:   default:
 72:     return jac->schur_user ? jac->schur_user : jac->pmat[1];
 73:   }
 74: }


 77: #include <petscdraw.h>
 80: static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
 81: {
 82:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
 83:   PetscErrorCode    ierr;
 84:   PetscBool         iascii,isdraw;
 85:   PetscInt          i,j;
 86:   PC_FieldSplitLink ilink = jac->head;

 89:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 90:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
 91:   if (iascii) {
 92:     if (jac->bs > 0) {
 93:       PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);
 94:     } else {
 95:       PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);
 96:     }
 97:     if (pc->useAmat) {
 98:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");
 99:     }
100:     if (jac->diag_use_amat) {
101:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for diagonal blocks\n");
102:     }
103:     if (jac->offdiag_use_amat) {
104:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for off-diagonal blocks\n");
105:     }
106:     PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");
107:     PetscViewerASCIIPushTab(viewer);
108:     for (i=0; i<jac->nsplits; i++) {
109:       if (ilink->fields) {
110:         PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
111:         PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
112:         for (j=0; j<ilink->nfields; j++) {
113:           if (j > 0) {
114:             PetscViewerASCIIPrintf(viewer,",");
115:           }
116:           PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
117:         }
118:         PetscViewerASCIIPrintf(viewer,"\n");
119:         PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
120:       } else {
121:         PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);
122:       }
123:       KSPView(ilink->ksp,viewer);
124:       ilink = ilink->next;
125:     }
126:     PetscViewerASCIIPopTab(viewer);
127:   }

129:  if (isdraw) {
130:     PetscDraw draw;
131:     PetscReal x,y,w,wd;

133:     PetscViewerDrawGetDraw(viewer,0,&draw);
134:     PetscDrawGetCurrentPoint(draw,&x,&y);
135:     w    = 2*PetscMin(1.0 - x,x);
136:     wd   = w/(jac->nsplits + 1);
137:     x    = x - wd*(jac->nsplits-1)/2.0;
138:     for (i=0; i<jac->nsplits; i++) {
139:       PetscDrawPushCurrentPoint(draw,x,y);
140:       KSPView(ilink->ksp,viewer);
141:       PetscDrawPopCurrentPoint(draw);
142:       x    += wd;
143:       ilink = ilink->next;
144:     }
145:   }
146:   return(0);
147: }

151: static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
152: {
153:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
154:   PetscErrorCode    ierr;
155:   PetscBool         iascii,isdraw;
156:   PetscInt          i,j;
157:   PC_FieldSplitLink ilink = jac->head;

160:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
161:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
162:   if (iascii) {
163:     if (jac->bs > 0) {
164:       PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);
165:     } else {
166:       PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);
167:     }
168:     if (pc->useAmat) {
169:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");
170:     }
171:     switch (jac->schurpre) {
172:     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
173:       PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");break;
174:     case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
175:       PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses (lumped, if requested) A00's diagonal's inverse\n");break;
176:     case PC_FIELDSPLIT_SCHUR_PRE_A11:
177:       PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");break;
178:     case PC_FIELDSPLIT_SCHUR_PRE_FULL:
179:       PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the exact Schur complement\n");break;
180:     case PC_FIELDSPLIT_SCHUR_PRE_USER:
181:       if (jac->schur_user) {
182:         PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");
183:       } else {
184:         PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");
185:       }
186:       break;
187:     default:
188:       SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
189:     }
190:     PetscViewerASCIIPrintf(viewer,"  Split info:\n");
191:     PetscViewerASCIIPushTab(viewer);
192:     for (i=0; i<jac->nsplits; i++) {
193:       if (ilink->fields) {
194:         PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
195:         PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
196:         for (j=0; j<ilink->nfields; j++) {
197:           if (j > 0) {
198:             PetscViewerASCIIPrintf(viewer,",");
199:           }
200:           PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
201:         }
202:         PetscViewerASCIIPrintf(viewer,"\n");
203:         PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
204:       } else {
205:         PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);
206:       }
207:       ilink = ilink->next;
208:     }
209:     PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");
210:     PetscViewerASCIIPushTab(viewer);
211:     if (jac->head) {
212:       KSPView(jac->head->ksp,viewer);
213:     } else  {PetscViewerASCIIPrintf(viewer,"  not yet available\n");}
214:     PetscViewerASCIIPopTab(viewer);
215:     if (jac->head && jac->kspupper != jac->head->ksp) {
216:       PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");
217:       PetscViewerASCIIPushTab(viewer);
218:       if (jac->kspupper) {KSPView(jac->kspupper,viewer);}
219:       else {PetscViewerASCIIPrintf(viewer,"  not yet available\n");}
220:       PetscViewerASCIIPopTab(viewer);
221:     }
222:     PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");
223:     PetscViewerASCIIPushTab(viewer);
224:     if (jac->kspschur) {
225:       KSPView(jac->kspschur,viewer);
226:     } else {
227:       PetscViewerASCIIPrintf(viewer,"  not yet available\n");
228:     }
229:     PetscViewerASCIIPopTab(viewer);
230:     PetscViewerASCIIPopTab(viewer);
231:   } else if (isdraw && jac->head) {
232:     PetscDraw draw;
233:     PetscReal x,y,w,wd,h;
234:     PetscInt  cnt = 2;
235:     char      str[32];

237:     PetscViewerDrawGetDraw(viewer,0,&draw);
238:     PetscDrawGetCurrentPoint(draw,&x,&y);
239:     if (jac->kspupper != jac->head->ksp) cnt++;
240:     w  = 2*PetscMin(1.0 - x,x);
241:     wd = w/(cnt + 1);

243:     PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);
244:     PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
245:     y   -= h;
246:     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
247:       PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);
248:     } else {
249:       PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);
250:     }
251:     PetscDrawStringBoxed(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
252:     y   -= h;
253:     x    = x - wd*(cnt-1)/2.0;

255:     PetscDrawPushCurrentPoint(draw,x,y);
256:     KSPView(jac->head->ksp,viewer);
257:     PetscDrawPopCurrentPoint(draw);
258:     if (jac->kspupper != jac->head->ksp) {
259:       x   += wd;
260:       PetscDrawPushCurrentPoint(draw,x,y);
261:       KSPView(jac->kspupper,viewer);
262:       PetscDrawPopCurrentPoint(draw);
263:     }
264:     x   += wd;
265:     PetscDrawPushCurrentPoint(draw,x,y);
266:     KSPView(jac->kspschur,viewer);
267:     PetscDrawPopCurrentPoint(draw);
268:   }
269:   return(0);
270: }

274: /* Precondition: jac->bs is set to a meaningful value */
275: static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
276: {
278:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
279:   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
280:   PetscBool      flg,flg_col;
281:   char           optionname[128],splitname[8],optionname_col[128];

284:   PetscMalloc1(jac->bs,&ifields);
285:   PetscMalloc1(jac->bs,&ifields_col);
286:   for (i=0,flg=PETSC_TRUE;; i++) {
287:     PetscSNPrintf(splitname,sizeof(splitname),"%D",i);
288:     PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);
289:     PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);
290:     nfields     = jac->bs;
291:     nfields_col = jac->bs;
292:     PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);
293:     PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);
294:     if (!flg) break;
295:     else if (flg && !flg_col) {
296:       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
297:       PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);
298:     } else {
299:       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
300:       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
301:       PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);
302:     }
303:   }
304:   if (i > 0) {
305:     /* Makes command-line setting of splits take precedence over setting them in code.
306:        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
307:        create new splits, which would probably not be what the user wanted. */
308:     jac->splitdefined = PETSC_TRUE;
309:   }
310:   PetscFree(ifields);
311:   PetscFree(ifields_col);
312:   return(0);
313: }

317: static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
318: {
319:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
320:   PetscErrorCode    ierr;
321:   PC_FieldSplitLink ilink = jac->head;
322:   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE,coupling = PETSC_FALSE;
323:   PetscInt          i;

326:   /*
327:    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
328:    Should probably be rewritten.
329:    */
330:   if (!ilink || jac->reset) {
331:     PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);
332:     PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_coupling",&coupling,NULL);
333:     if (pc->dm && jac->dm_splits && !stokes && !coupling) {
334:       PetscInt  numFields, f, i, j;
335:       char      **fieldNames;
336:       IS        *fields;
337:       DM        *dms;
338:       DM        subdm[128];
339:       PetscBool flg;

341:       DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);
342:       /* Allow the user to prescribe the splits */
343:       for (i = 0, flg = PETSC_TRUE;; i++) {
344:         PetscInt ifields[128];
345:         IS       compField;
346:         char     optionname[128], splitname[8];
347:         PetscInt nfields = numFields;

349:         PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);
350:         PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg);
351:         if (!flg) break;
352:         if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
353:         DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);
354:         if (nfields == 1) {
355:           PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);
356:         } else {
357:           PetscSNPrintf(splitname, sizeof(splitname), "%D", i);
358:           PCFieldSplitSetIS(pc, splitname, compField);
359:         }
360:         ISDestroy(&compField);
361:         for (j = 0; j < nfields; ++j) {
362:           f    = ifields[j];
363:           PetscFree(fieldNames[f]);
364:           ISDestroy(&fields[f]);
365:         }
366:       }
367:       if (i == 0) {
368:         for (f = 0; f < numFields; ++f) {
369:           PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);
370:           PetscFree(fieldNames[f]);
371:           ISDestroy(&fields[f]);
372:         }
373:       } else {
374:         for (j=0; j<numFields; j++) {
375:           DMDestroy(dms+j);
376:         }
377:         PetscFree(dms);
378:         PetscMalloc1(i, &dms);
379:         for (j = 0; j < i; ++j) dms[j] = subdm[j];
380:       }
381:       PetscFree(fieldNames);
382:       PetscFree(fields);
383:       if (dms) {
384:         PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");
385:         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
386:           const char *prefix;
387:           PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);
388:           PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);
389:           KSPSetDM(ilink->ksp, dms[i]);
390:           KSPSetDMActive(ilink->ksp, PETSC_FALSE);
391:           PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);
392:           DMDestroy(&dms[i]);
393:         }
394:         PetscFree(dms);
395:       }
396:     } else {
397:       if (jac->bs <= 0) {
398:         if (pc->pmat) {
399:           MatGetBlockSize(pc->pmat,&jac->bs);
400:         } else jac->bs = 1;
401:       }

403:       if (stokes) {
404:         IS       zerodiags,rest;
405:         PetscInt nmin,nmax;

407:         MatGetOwnershipRange(pc->mat,&nmin,&nmax);
408:         MatFindZeroDiagonals(pc->mat,&zerodiags);
409:         ISComplement(zerodiags,nmin,nmax,&rest);
410:         if (jac->reset) {
411:           if (!jac->head) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"At reset jac must have head from previous setup");
412:           jac->head->is       = rest;
413:           jac->head->next->is = zerodiags;
414:         } else {
415:           PCFieldSplitSetIS(pc,"0",rest);
416:           PCFieldSplitSetIS(pc,"1",zerodiags);
417:         }
418:         ISDestroy(&zerodiags);
419:         ISDestroy(&rest);
420:       } else if (coupling) {
421:         IS       coupling,rest;
422:         PetscInt nmin,nmax;

424:         MatGetOwnershipRange(pc->mat,&nmin,&nmax);
425:         MatFindOffBlockDiagonalEntries(pc->mat,&coupling);
426:         ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);
427:         ISSetIdentity(rest);
428:         if (jac->reset) {
429:           if (!jac->head) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"At reset jac must have head from previous setup");
430:           jac->head->is       = rest;
431:           jac->head->next->is = coupling;
432:         } else {
433:           PCFieldSplitSetIS(pc,"0",rest);
434:           PCFieldSplitSetIS(pc,"1",coupling);
435:         }
436:         ISDestroy(&coupling);
437:         ISDestroy(&rest);
438:       } else {
439:         if (jac->reset && !jac->isrestrict) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
440:         PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);
441:         if (!fieldsplit_default) {
442:           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
443:            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
444:           PCFieldSplitSetRuntimeSplits_Private(pc);
445:           if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
446:         }
447:         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
448:           PetscInfo(pc,"Using default splitting of fields\n");
449:           for (i=0; i<jac->bs; i++) {
450:             char splitname[8];
451:             PetscSNPrintf(splitname,sizeof(splitname),"%D",i);
452:             PCFieldSplitSetFields(pc,splitname,1,&i,&i);
453:           }
454:           jac->defaultsplit = PETSC_TRUE;
455:         }
456:       }
457:     }
458:   } else if (jac->nsplits == 1) {
459:     if (ilink->is) {
460:       IS       is2;
461:       PetscInt nmin,nmax;

463:       MatGetOwnershipRange(pc->mat,&nmin,&nmax);
464:       ISComplement(ilink->is,nmin,nmax,&is2);
465:       PCFieldSplitSetIS(pc,"1",is2);
466:       ISDestroy(&is2);
467:     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
468:   }

470:   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
471:   return(0);
472: }

474: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions,const char pre[], const char name[], char *value[], PetscBool *flg);

478: static PetscErrorCode PCSetUp_FieldSplit(PC pc)
479: {
480:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
481:   PetscErrorCode    ierr;
482:   PC_FieldSplitLink ilink;
483:   PetscInt          i,nsplit;
484:   PetscBool         sorted, sorted_col;

487:   PCFieldSplitSetDefaults(pc);
488:   nsplit = jac->nsplits;
489:   ilink  = jac->head;

491:   /* get the matrices for each split */
492:   if (!jac->issetup) {
493:     PetscInt rstart,rend,nslots,bs;

495:     jac->issetup = PETSC_TRUE;

497:     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
498:     if (jac->defaultsplit || !ilink->is) {
499:       if (jac->bs <= 0) jac->bs = nsplit;
500:     }
501:     bs     = jac->bs;
502:     MatGetOwnershipRange(pc->pmat,&rstart,&rend);
503:     nslots = (rend - rstart)/bs;
504:     for (i=0; i<nsplit; i++) {
505:       if (jac->defaultsplit) {
506:         ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);
507:         ISDuplicate(ilink->is,&ilink->is_col);
508:       } else if (!ilink->is) {
509:         if (ilink->nfields > 1) {
510:           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
511:           PetscMalloc1(ilink->nfields*nslots,&ii);
512:           PetscMalloc1(ilink->nfields*nslots,&jj);
513:           for (j=0; j<nslots; j++) {
514:             for (k=0; k<nfields; k++) {
515:               ii[nfields*j + k] = rstart + bs*j + fields[k];
516:               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
517:             }
518:           }
519:           ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);
520:           ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);
521:           ISSetBlockSize(ilink->is, nfields);
522:           ISSetBlockSize(ilink->is_col, nfields);
523:         } else {
524:           ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);
525:           ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);
526:        }
527:       }
528:       ISSorted(ilink->is,&sorted);
529:       if (ilink->is_col) { ISSorted(ilink->is_col,&sorted_col); }
530:       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
531:       ilink = ilink->next;
532:     }
533:   }

535:   ilink = jac->head;
536:   if (!jac->pmat) {
537:     Vec xtmp;

539:     MatCreateVecs(pc->pmat,&xtmp,NULL);
540:     PetscMalloc1(nsplit,&jac->pmat);
541:     PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);
542:     for (i=0; i<nsplit; i++) {
543:       MatNullSpace sp;

545:       /* Check for preconditioning matrix attached to IS */
546:       PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);
547:       if (jac->pmat[i]) {
548:         PetscObjectReference((PetscObject) jac->pmat[i]);
549:         if (jac->type == PC_COMPOSITE_SCHUR) {
550:           jac->schur_user = jac->pmat[i];

552:           PetscObjectReference((PetscObject) jac->schur_user);
553:         }
554:       } else {
555:         const char *prefix;
556:         MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);
557:         KSPGetOptionsPrefix(ilink->ksp,&prefix);
558:         MatSetOptionsPrefix(jac->pmat[i],prefix);
559:         MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");
560:       }
561:       /* create work vectors for each split */
562:       MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);
563:       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
564:       /* compute scatter contexts needed by multiplicative versions and non-default splits */
565:       VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);
566:       PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);
567:       if (sp) {
568:         MatSetNearNullSpace(jac->pmat[i], sp);
569:       }
570:       ilink = ilink->next;
571:     }
572:     VecDestroy(&xtmp);
573:   } else {
574:     for (i=0; i<nsplit; i++) {
575:       Mat pmat;

577:       /* Check for preconditioning matrix attached to IS */
578:       PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);
579:       if (!pmat) {
580:         MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);
581:       }
582:       ilink = ilink->next;
583:     }
584:   }
585:   if (jac->diag_use_amat) {
586:     ilink = jac->head;
587:     if (!jac->mat) {
588:       PetscMalloc1(nsplit,&jac->mat);
589:       for (i=0; i<nsplit; i++) {
590:         MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);
591:         ilink = ilink->next;
592:       }
593:     } else {
594:       for (i=0; i<nsplit; i++) {
595:         if (jac->mat[i]) {MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);}
596:         ilink = ilink->next;
597:       }
598:     }
599:   } else {
600:     jac->mat = jac->pmat;
601:   }

603:   /* Check for null space attached to IS */
604:   ilink = jac->head;
605:   for (i=0; i<nsplit; i++) {
606:     MatNullSpace sp;

608:     PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);
609:     if (sp) {
610:       MatSetNullSpace(jac->mat[i], sp);
611:     }
612:     ilink = ilink->next;
613:   }

615:   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
616:     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
617:     /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
618:     ilink = jac->head;
619:     if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
620:       /* special case need where Afield[0] is not needed and only certain columns of Afield[1] are needed since update is only on those rows of the solution */
621:       if (!jac->Afield) {
622:         PetscCalloc1(nsplit,&jac->Afield);
623:         MatGetSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);
624:       } else {
625:         MatGetSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_REUSE_MATRIX,&jac->Afield[1]);
626:       }
627:     } else {
628:       if (!jac->Afield) {
629:         PetscMalloc1(nsplit,&jac->Afield);
630:         for (i=0; i<nsplit; i++) {
631:           MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);
632:           ilink = ilink->next;
633:         }
634:       } else {
635:         for (i=0; i<nsplit; i++) {
636:           MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);
637:           ilink = ilink->next;
638:         }
639:       }
640:     }
641:   }

643:   if (jac->type == PC_COMPOSITE_SCHUR) {
644:     IS          ccis;
645:     PetscInt    rstart,rend;
646:     char        lscname[256];
647:     PetscObject LSC_L;

649:     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");

651:     /* When extracting off-diagonal submatrices, we take complements from this range */
652:     MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);

654:     /* need to handle case when one is resetting up the preconditioner */
655:     if (jac->schur) {
656:       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;

658:       MatSchurComplementGetKSP(jac->schur, &kspInner);
659:       ilink = jac->head;
660:       ISComplement(ilink->is_col,rstart,rend,&ccis);
661:       if (jac->offdiag_use_amat) {
662:         MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);
663:       } else {
664:         MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);
665:       }
666:       ISDestroy(&ccis);
667:       ilink = ilink->next;
668:       ISComplement(ilink->is_col,rstart,rend,&ccis);
669:       if (jac->offdiag_use_amat) {
670:         MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);
671:       } else {
672:         MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);
673:       }
674:       ISDestroy(&ccis);
675:       MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
676:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
677:         MatDestroy(&jac->schurp);
678:         MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);
679:       }
680:       if (kspA != kspInner) {
681:         KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);
682:       }
683:       if (kspUpper != kspA) {
684:         KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);
685:       }
686:       KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));
687:     } else {
688:       const char   *Dprefix;
689:       char         schurprefix[256], schurmatprefix[256];
690:       char         schurtestoption[256];
691:       MatNullSpace sp;
692:       PetscBool    flg;

694:       /* extract the A01 and A10 matrices */
695:       ilink = jac->head;
696:       ISComplement(ilink->is_col,rstart,rend,&ccis);
697:       if (jac->offdiag_use_amat) {
698:         MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
699:       } else {
700:         MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
701:       }
702:       ISDestroy(&ccis);
703:       ilink = ilink->next;
704:       ISComplement(ilink->is_col,rstart,rend,&ccis);
705:       if (jac->offdiag_use_amat) {
706:         MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
707:       } else {
708:         MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
709:       }
710:       ISDestroy(&ccis);

712:       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
713:       MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);
714:       MatSetType(jac->schur,MATSCHURCOMPLEMENT);
715:       MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
716:       PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
717:       /* Note that the inner KSP is NOT going to inherit this prefix, and if it did, it would be reset just below.  Is that what we want? */
718:       MatSetOptionsPrefix(jac->schur,schurmatprefix);
719:       MatSetFromOptions(jac->schur);
720:       MatGetNullSpace(jac->mat[1], &sp);
721:       if (sp) {
722:         MatSetNullSpace(jac->schur, sp);
723:       }

725:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);
726:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
727:       if (flg) {
728:         DM  dmInner;
729:         KSP kspInner;

731:         MatSchurComplementGetKSP(jac->schur, &kspInner);
732:         PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
733:         /* Indent this deeper to emphasize the "inner" nature of this solver. */
734:         PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);
735:         KSPSetOptionsPrefix(kspInner, schurprefix);

737:         /* Set DM for new solver */
738:         KSPGetDM(jac->head->ksp, &dmInner);
739:         KSPSetDM(kspInner, dmInner);
740:         KSPSetDMActive(kspInner, PETSC_FALSE);
741:       } else {
742:          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
743:           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
744:           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
745:           * S = D - C A_inner^{-1} B, we expect S to be defined using an accurate definition of A_inner^{-1}, so we make
746:           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
747:           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
748:         KSPSetType(jac->head->ksp,KSPGMRES);
749:         MatSchurComplementSetKSP(jac->schur,jac->head->ksp);
750:       }
751:       KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);
752:       KSPSetFromOptions(jac->head->ksp);
753:       MatSetFromOptions(jac->schur);

755:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);
756:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
757:       if (flg) {
758:         DM dmInner;

760:         PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
761:         KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);
762:         KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);
763:         KSPSetOptionsPrefix(jac->kspupper, schurprefix);
764:         KSPGetDM(jac->head->ksp, &dmInner);
765:         KSPSetDM(jac->kspupper, dmInner);
766:         KSPSetDMActive(jac->kspupper, PETSC_FALSE);
767:         KSPSetFromOptions(jac->kspupper);
768:         KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);
769:         VecDuplicate(jac->head->x, &jac->head->z);
770:       } else {
771:         jac->kspupper = jac->head->ksp;
772:         PetscObjectReference((PetscObject) jac->head->ksp);
773:       }

775:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
776:         MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);
777:       }
778:       KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);
779:       KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);
780:       PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);
781:       PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);
782:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
783:         PC pcschur;
784:         KSPGetPC(jac->kspschur,&pcschur);
785:         PCSetType(pcschur,PCNONE);
786:         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
787:       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
788:         MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);
789:       }
790:       KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));
791:       KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);
792:       KSPSetOptionsPrefix(jac->kspschur,         Dprefix);
793:       /* propogate DM */
794:       {
795:         DM sdm;
796:         KSPGetDM(jac->head->next->ksp, &sdm);
797:         if (sdm) {
798:           KSPSetDM(jac->kspschur, sdm);
799:           KSPSetDMActive(jac->kspschur, PETSC_FALSE);
800:         }
801:       }
802:       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
803:       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
804:       KSPSetFromOptions(jac->kspschur);
805:     }

807:     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
808:     PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);
809:     PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);
810:     if (!LSC_L) {PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);}
811:     if (LSC_L) {PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);}
812:     PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);
813:     PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);
814:     if (!LSC_L) {PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);}
815:     if (LSC_L) {PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);}
816:   } else {
817:     /* set up the individual splits' PCs */
818:     i     = 0;
819:     ilink = jac->head;
820:     while (ilink) {
821:       KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);
822:       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
823:       if (!jac->suboptionsset) {KSPSetFromOptions(ilink->ksp);}
824:       i++;
825:       ilink = ilink->next;
826:     }
827:   }

829:   jac->suboptionsset = PETSC_TRUE;
830:   return(0);
831: }

833: #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
834:   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
835:    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
836:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
837:    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
838:    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
839:    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
840:    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))

844: static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
845: {
846:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
847:   PetscErrorCode    ierr;
848:   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
849:   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;

852:   switch (jac->schurfactorization) {
853:   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
854:     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
855:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
856:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
857:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
858:     PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
859:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
860:     PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
861:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
862:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
863:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
864:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
865:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
866:     VecScale(ilinkD->y,-1.);
867:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
868:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
869:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
870:     break;
871:   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
872:     /* [A00 0; A10 S], suitable for left preconditioning */
873:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
874:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
875:     PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
876:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
877:     PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
878:     MatMult(jac->C,ilinkA->y,ilinkD->x);
879:     VecScale(ilinkD->x,-1.);
880:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
881:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
882:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
883:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
884:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
885:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
886:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
887:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
888:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
889:     break;
890:   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
891:     /* [A00 A01; 0 S], suitable for right preconditioning */
892:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
893:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
894:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
895:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
896:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);    MatMult(jac->B,ilinkD->y,ilinkA->x);
897:     VecScale(ilinkA->x,-1.);
898:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
899:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
900:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
901:     PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
902:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
903:     PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
904:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
905:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
906:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
907:     break;
908:   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
909:     /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1], an exact solve if applied exactly, needs one extra solve with A */
910:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
911:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
912:     PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);
913:     KSPSolve(kspLower,ilinkA->x,ilinkA->y);
914:     PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);
915:     MatMult(jac->C,ilinkA->y,ilinkD->x);
916:     VecScale(ilinkD->x,-1.0);
917:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
918:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);

920:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
921:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
922:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
923:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
924:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);

926:     if (kspUpper == kspA) {
927:       MatMult(jac->B,ilinkD->y,ilinkA->y);
928:       VecAXPY(ilinkA->x,-1.0,ilinkA->y);
929:       PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
930:       KSPSolve(kspA,ilinkA->x,ilinkA->y);
931:       PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
932:     } else {
933:       PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
934:       KSPSolve(kspA,ilinkA->x,ilinkA->y);
935:       PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
936:       MatMult(jac->B,ilinkD->y,ilinkA->x);
937:       PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);
938:       KSPSolve(kspUpper,ilinkA->x,ilinkA->z);
939:       PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);
940:       VecAXPY(ilinkA->y,-1.0,ilinkA->z);
941:     }
942:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
943:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
944:   }
945:   return(0);
946: }

950: static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
951: {
952:   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
953:   PetscErrorCode     ierr;
954:   PC_FieldSplitLink  ilink = jac->head;
955:   PetscInt           cnt,bs;
956:   KSPConvergedReason reason;

959:   if (jac->type == PC_COMPOSITE_ADDITIVE) {
960:     if (jac->defaultsplit) {
961:       VecGetBlockSize(x,&bs);
962:       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
963:       VecGetBlockSize(y,&bs);
964:       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
965:       VecStrideGatherAll(x,jac->x,INSERT_VALUES);
966:       while (ilink) {
967:         PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
968:         KSPSolve(ilink->ksp,ilink->x,ilink->y);
969:         PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
970:         KSPGetConvergedReason(ilink->ksp,&reason);
971:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
972:           pc->failedreason = PC_SUBPC_ERROR;
973:         }
974:         ilink = ilink->next;
975:       }
976:       VecStrideScatterAll(jac->y,y,INSERT_VALUES);
977:     } else {
978:       VecSet(y,0.0);
979:       while (ilink) {
980:         FieldSplitSplitSolveAdd(ilink,x,y);
981:         KSPGetConvergedReason(ilink->ksp,&reason);
982:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
983:           pc->failedreason = PC_SUBPC_ERROR;
984:         }
985:         ilink = ilink->next;
986:       }
987:     }
988:   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
989:     VecSet(y,0.0);
990:     /* solve on first block for first block variables */
991:     VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
992:     VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
993:     PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
994:     KSPSolve(ilink->ksp,ilink->x,ilink->y);
995:     PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
996:     KSPGetConvergedReason(ilink->ksp,&reason);
997:     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
998:       pc->failedreason = PC_SUBPC_ERROR;
999:     }
1000:     VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1001:     VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);

1003:     /* compute the residual only onto second block variables using first block variables */
1004:     MatMult(jac->Afield[1],ilink->y,ilink->next->x);
1005:     ilink = ilink->next;
1006:     VecScale(ilink->x,-1.0);
1007:     VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1008:     VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);

1010:     /* solve on second block variables */
1011:     PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1012:     KSPSolve(ilink->ksp,ilink->x,ilink->y);
1013:     PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1014:     KSPGetConvergedReason(ilink->ksp,&reason);
1015:     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1016:       pc->failedreason = PC_SUBPC_ERROR;
1017:     }
1018:     VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1019:     VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1020:   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1021:     if (!jac->w1) {
1022:       VecDuplicate(x,&jac->w1);
1023:       VecDuplicate(x,&jac->w2);
1024:     }
1025:     VecSet(y,0.0);
1026:     FieldSplitSplitSolveAdd(ilink,x,y);
1027:     KSPGetConvergedReason(ilink->ksp,&reason);
1028:     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1029:       pc->failedreason = PC_SUBPC_ERROR;
1030:     }
1031:     cnt  = 1;
1032:     while (ilink->next) {
1033:       ilink = ilink->next;
1034:       /* compute the residual only over the part of the vector needed */
1035:       MatMult(jac->Afield[cnt++],y,ilink->x);
1036:       VecScale(ilink->x,-1.0);
1037:       VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1038:       VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1039:       PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1040:       KSPSolve(ilink->ksp,ilink->x,ilink->y);
1041:       PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1042:       KSPGetConvergedReason(ilink->ksp,&reason);
1043:       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1044:         pc->failedreason = PC_SUBPC_ERROR;
1045:       }
1046:       VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1047:       VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1048:     }
1049:     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1050:       cnt -= 2;
1051:       while (ilink->previous) {
1052:         ilink = ilink->previous;
1053:         /* compute the residual only over the part of the vector needed */
1054:         MatMult(jac->Afield[cnt--],y,ilink->x);
1055:         VecScale(ilink->x,-1.0);
1056:         VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1057:         VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1058:         PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1059:         KSPSolve(ilink->ksp,ilink->x,ilink->y);
1060:         PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1061:         KSPGetConvergedReason(ilink->ksp,&reason);
1062:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1063:           pc->failedreason = PC_SUBPC_ERROR;
1064:         }
1065:         VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1066:         VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1067:       }
1068:     }
1069:   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
1070:   return(0);
1071: }

1073: #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1074:   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1075:    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1076:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1077:    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) ||                  \
1078:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1079:    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1080:    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))

1084: static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1085: {
1086:   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1087:   PetscErrorCode     ierr;
1088:   PC_FieldSplitLink  ilink = jac->head;
1089:   PetscInt           bs;
1090:   KSPConvergedReason reason;

1093:   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1094:     if (jac->defaultsplit) {
1095:       VecGetBlockSize(x,&bs);
1096:       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1097:       VecGetBlockSize(y,&bs);
1098:       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1099:       VecStrideGatherAll(x,jac->x,INSERT_VALUES);
1100:       while (ilink) {
1101:         PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1102:         KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);
1103:         PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1104:         KSPGetConvergedReason(ilink->ksp,&reason);
1105:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1106:           pc->failedreason = PC_SUBPC_ERROR;
1107:         }
1108:         ilink = ilink->next;
1109:       }
1110:       VecStrideScatterAll(jac->y,y,INSERT_VALUES);
1111:     } else {
1112:       VecSet(y,0.0);
1113:       while (ilink) {
1114:         FieldSplitSplitSolveAddTranspose(ilink,x,y);
1115:         KSPGetConvergedReason(ilink->ksp,&reason);
1116:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1117:           pc->failedreason = PC_SUBPC_ERROR;
1118:         }
1119:         ilink = ilink->next;
1120:       }
1121:     }
1122:   } else {
1123:     if (!jac->w1) {
1124:       VecDuplicate(x,&jac->w1);
1125:       VecDuplicate(x,&jac->w2);
1126:     }
1127:     VecSet(y,0.0);
1128:     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1129:       FieldSplitSplitSolveAddTranspose(ilink,x,y);
1130:       KSPGetConvergedReason(ilink->ksp,&reason);
1131:       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1132:         pc->failedreason = PC_SUBPC_ERROR;
1133:       }
1134:       while (ilink->next) {
1135:         ilink = ilink->next;
1136:         MatMultTranspose(pc->mat,y,jac->w1);
1137:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1138:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1139:       }
1140:       while (ilink->previous) {
1141:         ilink = ilink->previous;
1142:         MatMultTranspose(pc->mat,y,jac->w1);
1143:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1144:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1145:       }
1146:     } else {
1147:       while (ilink->next) {   /* get to last entry in linked list */
1148:         ilink = ilink->next;
1149:       }
1150:       FieldSplitSplitSolveAddTranspose(ilink,x,y);
1151:       KSPGetConvergedReason(ilink->ksp,&reason);
1152:       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1153:         pc->failedreason = PC_SUBPC_ERROR;
1154:       }
1155:       while (ilink->previous) {
1156:         ilink = ilink->previous;
1157:         MatMultTranspose(pc->mat,y,jac->w1);
1158:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1159:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1160:       }
1161:     }
1162:   }
1163:   return(0);
1164: }

1168: static PetscErrorCode PCReset_FieldSplit(PC pc)
1169: {
1170:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1171:   PetscErrorCode    ierr;
1172:   PC_FieldSplitLink ilink = jac->head,next;

1175:   while (ilink) {
1176:     KSPReset(ilink->ksp);
1177:     VecDestroy(&ilink->x);
1178:     VecDestroy(&ilink->y);
1179:     VecDestroy(&ilink->z);
1180:     VecScatterDestroy(&ilink->sctx);
1181:     if (!ilink->is_orig) {             /* save the original IS */
1182:       PetscObjectReference((PetscObject)ilink->is);
1183:       ilink->is_orig = ilink->is;
1184:     }
1185:     ISDestroy(&ilink->is);
1186:     ISDestroy(&ilink->is_col);
1187:     next  = ilink->next;
1188:     ilink = next;
1189:   }
1190:   PetscFree2(jac->x,jac->y);
1191:   if (jac->mat && jac->mat != jac->pmat) {
1192:     MatDestroyMatrices(jac->nsplits,&jac->mat);
1193:   } else if (jac->mat) {
1194:     jac->mat = NULL;
1195:   }
1196:   if (jac->pmat) {MatDestroyMatrices(jac->nsplits,&jac->pmat);}
1197:   if (jac->Afield) {MatDestroyMatrices(jac->nsplits,&jac->Afield);}
1198:   VecDestroy(&jac->w1);
1199:   VecDestroy(&jac->w2);
1200:   MatDestroy(&jac->schur);
1201:   MatDestroy(&jac->schurp);
1202:   MatDestroy(&jac->schur_user);
1203:   KSPDestroy(&jac->kspschur);
1204:   KSPDestroy(&jac->kspupper);
1205:   MatDestroy(&jac->B);
1206:   MatDestroy(&jac->C);
1207:   jac->reset = PETSC_TRUE;
1208:   jac->isrestrict = PETSC_FALSE;
1209:   return(0);
1210: }

1214: static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1215: {
1216:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1217:   PetscErrorCode    ierr;
1218:   PC_FieldSplitLink ilink = jac->head,next;

1221:   PCReset_FieldSplit(pc);
1222:   while (ilink) {
1223:     KSPDestroy(&ilink->ksp);
1224:     ISDestroy(&ilink->is_orig);
1225:     next  = ilink->next;
1226:     PetscFree(ilink->splitname);
1227:     PetscFree(ilink->fields);
1228:     PetscFree(ilink->fields_col);
1229:     PetscFree(ilink);
1230:     ilink = next;
1231:   }
1232:   PetscFree2(jac->x,jac->y);
1233:   PetscFree(pc->data);
1234:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);
1235:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);
1236:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);
1237:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);
1238:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);
1239:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);
1240:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);
1241:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);
1242:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);
1243:   return(0);
1244: }

1248: static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
1249: {
1250:   PetscErrorCode  ierr;
1251:   PetscInt        bs;
1252:   PetscBool       flg,stokes = PETSC_FALSE;
1253:   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
1254:   PCCompositeType ctype;

1257:   PetscOptionsHead(PetscOptionsObject,"FieldSplit options");
1258:   PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);
1259:   PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);
1260:   if (flg) {
1261:     PCFieldSplitSetBlockSize(pc,bs);
1262:   }
1263:   jac->diag_use_amat = pc->useAmat;
1264:   PetscOptionsBool("-pc_fieldsplit_diag_use_amat","Use Amat (not Pmat) to extract diagonal fieldsplit blocks", "PCFieldSplitSetDiagUseAmat",jac->diag_use_amat,&jac->diag_use_amat,NULL);
1265:   jac->offdiag_use_amat = pc->useAmat;
1266:   PetscOptionsBool("-pc_fieldsplit_off_diag_use_amat","Use Amat (not Pmat) to extract off-diagonal fieldsplit blocks", "PCFieldSplitSetOffDiagUseAmat",jac->offdiag_use_amat,&jac->offdiag_use_amat,NULL);
1267:   /* FIXME: No programmatic equivalent to the following. */
1268:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);
1269:   if (stokes) {
1270:     PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
1271:     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1272:   }

1274:   PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);
1275:   if (flg) {
1276:     PCFieldSplitSetType(pc,ctype);
1277:   }
1278:   /* Only setup fields once */
1279:   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1280:     /* only allow user to set fields from command line if bs is already known.
1281:        otherwise user can set them in PCFieldSplitSetDefaults() */
1282:     PCFieldSplitSetRuntimeSplits_Private(pc);
1283:     if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
1284:   }
1285:   if (jac->type == PC_COMPOSITE_SCHUR) {
1286:     PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);
1287:     if (flg) {PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");}
1288:     PetscOptionsEnum("-pc_fieldsplit_schur_fact_type","Which off-diagonal parts of the block factorization to use","PCFieldSplitSetSchurFactType",PCFieldSplitSchurFactTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,NULL);
1289:     PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);
1290:   }
1291:   PetscOptionsTail();
1292:   return(0);
1293: }

1295: /*------------------------------------------------------------------------------------*/

1299: static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1300: {
1301:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1302:   PetscErrorCode    ierr;
1303:   PC_FieldSplitLink ilink,next = jac->head;
1304:   char              prefix[128];
1305:   PetscInt          i;

1308:   if (jac->splitdefined) {
1309:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1310:     return(0);
1311:   }
1312:   for (i=0; i<n; i++) {
1313:     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1314:     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1315:   }
1316:   PetscNew(&ilink);
1317:   if (splitname) {
1318:     PetscStrallocpy(splitname,&ilink->splitname);
1319:   } else {
1320:     PetscMalloc1(3,&ilink->splitname);
1321:     PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);
1322:   }
1323:   ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + jac->nsplits : KSP_Solve_FS_0 + 4; /* Any split great than 4 gets logged in the 4th split */
1324:   PetscMalloc1(n,&ilink->fields);
1325:   PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));
1326:   PetscMalloc1(n,&ilink->fields_col);
1327:   PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));

1329:   ilink->nfields = n;
1330:   ilink->next    = NULL;
1331:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1332:   KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1333:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1334:   KSPSetType(ilink->ksp,KSPPREONLY);
1335:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

1337:   PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);
1338:   KSPSetOptionsPrefix(ilink->ksp,prefix);

1340:   if (!next) {
1341:     jac->head       = ilink;
1342:     ilink->previous = NULL;
1343:   } else {
1344:     while (next->next) {
1345:       next = next->next;
1346:     }
1347:     next->next      = ilink;
1348:     ilink->previous = next;
1349:   }
1350:   jac->nsplits++;
1351:   return(0);
1352: }

1356: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1357: {
1358:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1362:   PetscMalloc1(jac->nsplits,subksp);
1363:   MatSchurComplementGetKSP(jac->schur,*subksp);

1365:   (*subksp)[1] = jac->kspschur;
1366:   if (n) *n = jac->nsplits;
1367:   return(0);
1368: }

1372: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1373: {
1374:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1375:   PetscErrorCode    ierr;
1376:   PetscInt          cnt   = 0;
1377:   PC_FieldSplitLink ilink = jac->head;

1380:   PetscMalloc1(jac->nsplits,subksp);
1381:   while (ilink) {
1382:     (*subksp)[cnt++] = ilink->ksp;
1383:     ilink            = ilink->next;
1384:   }
1385:   if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number of splits in linked list %D does not match number in object %D",cnt,jac->nsplits);
1386:   if (n) *n = jac->nsplits;
1387:   return(0);
1388: }

1392: /*@C
1393:     PCFieldSplitRestrictIS - Restricts the fieldsplit ISs to be within a given IS.

1395:     Input Parameters:
1396: +   pc  - the preconditioner context
1397: +   is - the index set that defines the indices to which the fieldsplit is to be restricted

1399:     Level: advanced

1401: @*/
1402: PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
1403: {

1409:   PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));
1410:   return(0);
1411: }


1416: static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
1417: {
1418:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1419:   PetscErrorCode    ierr;
1420:   PC_FieldSplitLink ilink = jac->head, next;
1421:   PetscInt          localsize,size,sizez,i;
1422:   const PetscInt    *ind, *indz;
1423:   PetscInt          *indc, *indcz;
1424:   PetscBool         flg;

1427:   ISGetLocalSize(isy,&localsize);
1428:   MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));
1429:   size -= localsize;
1430:   while(ilink) {
1431:     IS isrl,isr;
1432:     PC subpc;
1433:     if (jac->reset) {
1434:       ISEmbed(ilink->is_orig, isy, PETSC_TRUE, &isrl);
1435:     } else {
1436:       ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);
1437:     }
1438:     ISGetLocalSize(isrl,&localsize);
1439:     PetscMalloc1(localsize,&indc);
1440:     ISGetIndices(isrl,&ind);
1441:     PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));
1442:     ISRestoreIndices(isrl,&ind);
1443:     ISDestroy(&isrl);
1444:     for (i=0; i<localsize; i++) *(indc+i) += size;
1445:     ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);
1446:     PetscObjectReference((PetscObject)isr);
1447:     ISDestroy(&ilink->is);
1448:     ilink->is     = isr;
1449:     PetscObjectReference((PetscObject)isr);
1450:     ISDestroy(&ilink->is_col);
1451:     ilink->is_col = isr;
1452:     ISDestroy(&isr);
1453:     KSPGetPC(ilink->ksp, &subpc);
1454:     PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);
1455:     if(flg) {
1456:       IS iszl,isz;
1457:       MPI_Comm comm;
1458:       if (jac->reset) {
1459:         ISGetLocalSize(ilink->is_orig,&localsize);
1460:         comm = PetscObjectComm((PetscObject)ilink->is_orig);
1461:         ISEmbed(isy, ilink->is_orig, PETSC_TRUE, &iszl);
1462:       } else {
1463:         ISGetLocalSize(ilink->is,&localsize);
1464:         comm = PetscObjectComm((PetscObject)ilink->is);
1465:         ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);
1466:       }
1467:       MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);
1468:       sizez -= localsize;
1469:       ISGetLocalSize(iszl,&localsize);
1470:       PetscMalloc1(localsize,&indcz);
1471:       ISGetIndices(iszl,&indz);
1472:       PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));
1473:       ISRestoreIndices(iszl,&indz);
1474:       ISDestroy(&iszl);
1475:       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
1476:       ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);
1477:       PCFieldSplitRestrictIS(subpc,isz);
1478:       ISDestroy(&isz);
1479:     }
1480:     next = ilink->next;
1481:     ilink = next;
1482:   }
1483:   jac->isrestrict = PETSC_TRUE;
1484:   return(0);
1485: }

1489: static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1490: {
1491:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1492:   PetscErrorCode    ierr;
1493:   PC_FieldSplitLink ilink, next = jac->head;
1494:   char              prefix[128];

1497:   if (jac->splitdefined) {
1498:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1499:     return(0);
1500:   }
1501:   PetscNew(&ilink);
1502:   if (splitname) {
1503:     PetscStrallocpy(splitname,&ilink->splitname);
1504:   } else {
1505:     PetscMalloc1(8,&ilink->splitname);
1506:     PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);
1507:   }
1508:   ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + jac->nsplits : KSP_Solve_FS_0 + 4; /* Any split great than 4 gets logged in the 4th split */
1509:   PetscObjectReference((PetscObject)is);
1510:   ISDestroy(&ilink->is);
1511:   ilink->is     = is;
1512:   PetscObjectReference((PetscObject)is);
1513:   ISDestroy(&ilink->is_col);
1514:   ilink->is_col = is;
1515:   ilink->next   = NULL;
1516:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1517:   KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1518:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1519:   KSPSetType(ilink->ksp,KSPPREONLY);
1520:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

1522:   PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);
1523:   KSPSetOptionsPrefix(ilink->ksp,prefix);

1525:   if (!next) {
1526:     jac->head       = ilink;
1527:     ilink->previous = NULL;
1528:   } else {
1529:     while (next->next) {
1530:       next = next->next;
1531:     }
1532:     next->next      = ilink;
1533:     ilink->previous = next;
1534:   }
1535:   jac->nsplits++;
1536:   return(0);
1537: }

1541: /*@
1542:     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner

1544:     Logically Collective on PC

1546:     Input Parameters:
1547: +   pc  - the preconditioner context
1548: .   splitname - name of this split, if NULL the number of the split is used
1549: .   n - the number of fields in this split
1550: -   fields - the fields in this split

1552:     Level: intermediate

1554:     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.

1556:      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1557:      size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean
1558:      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1559:      where the numbered entries indicate what is in the field.

1561:      This function is called once per split (it creates a new split each time).  Solve options
1562:      for this split will be available under the prefix -fieldsplit_SPLITNAME_.

1564:      Developer Note: This routine does not actually create the IS representing the split, that is delayed
1565:      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
1566:      available when this routine is called.

1568: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()

1570: @*/
1571: PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1572: {

1578:   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1580:   PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));
1581:   return(0);
1582: }

1586: /*@
1587:     PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)

1589:     Logically Collective on PC

1591:     Input Parameters:
1592: +   pc  - the preconditioner object
1593: -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from

1595:     Options Database:
1596: .     -pc_fieldsplit_diag_use_amat

1598:     Level: intermediate

1600: .seealso: PCFieldSplitGetDiagUseAmat(), PCFieldSplitSetOffDiagUseAmat(), PCFIELDSPLIT

1602: @*/
1603: PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1604: {
1605:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1606:   PetscBool      isfs;

1611:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1612:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1613:   jac->diag_use_amat = flg;
1614:   return(0);
1615: }

1619: /*@
1620:     PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)

1622:     Logically Collective on PC

1624:     Input Parameters:
1625: .   pc  - the preconditioner object

1627:     Output Parameters:
1628: .   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from


1631:     Level: intermediate

1633: .seealso: PCFieldSplitSetDiagUseAmat(), PCFieldSplitGetOffDiagUseAmat(), PCFIELDSPLIT

1635: @*/
1636: PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1637: {
1638:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1639:   PetscBool      isfs;

1645:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1646:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1647:   *flg = jac->diag_use_amat;
1648:   return(0);
1649: }

1653: /*@
1654:     PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)

1656:     Logically Collective on PC

1658:     Input Parameters:
1659: +   pc  - the preconditioner object
1660: -   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from

1662:     Options Database:
1663: .     -pc_fieldsplit_off_diag_use_amat

1665:     Level: intermediate

1667: .seealso: PCFieldSplitGetOffDiagUseAmat(), PCFieldSplitSetDiagUseAmat(), PCFIELDSPLIT

1669: @*/
1670: PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1671: {
1672:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1673:   PetscBool      isfs;

1678:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1679:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1680:   jac->offdiag_use_amat = flg;
1681:   return(0);
1682: }

1686: /*@
1687:     PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)

1689:     Logically Collective on PC

1691:     Input Parameters:
1692: .   pc  - the preconditioner object

1694:     Output Parameters:
1695: .   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from


1698:     Level: intermediate

1700: .seealso: PCFieldSplitSetOffDiagUseAmat(), PCFieldSplitGetDiagUseAmat(), PCFIELDSPLIT

1702: @*/
1703: PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1704: {
1705:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1706:   PetscBool      isfs;

1712:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1713:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1714:   *flg = jac->offdiag_use_amat;
1715:   return(0);
1716: }



1722: /*@C
1723:     PCFieldSplitSetIS - Sets the exact elements for field

1725:     Logically Collective on PC

1727:     Input Parameters:
1728: +   pc  - the preconditioner context
1729: .   splitname - name of this split, if NULL the number of the split is used
1730: -   is - the index set that defines the vector elements in this field


1733:     Notes:
1734:     Use PCFieldSplitSetFields(), for fields defined by strided types.

1736:     This function is called once per split (it creates a new split each time).  Solve options
1737:     for this split will be available under the prefix -fieldsplit_SPLITNAME_.

1739:     Level: intermediate

1741: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()

1743: @*/
1744: PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1745: {

1752:   PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));
1753:   return(0);
1754: }

1758: /*@
1759:     PCFieldSplitGetIS - Retrieves the elements for a field as an IS

1761:     Logically Collective on PC

1763:     Input Parameters:
1764: +   pc  - the preconditioner context
1765: -   splitname - name of this split

1767:     Output Parameter:
1768: -   is - the index set that defines the vector elements in this field, or NULL if the field is not found

1770:     Level: intermediate

1772: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()

1774: @*/
1775: PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
1776: {

1783:   {
1784:     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
1785:     PC_FieldSplitLink ilink = jac->head;
1786:     PetscBool         found;

1788:     *is = NULL;
1789:     while (ilink) {
1790:       PetscStrcmp(ilink->splitname, splitname, &found);
1791:       if (found) {
1792:         *is = ilink->is;
1793:         break;
1794:       }
1795:       ilink = ilink->next;
1796:     }
1797:   }
1798:   return(0);
1799: }

1803: /*@
1804:     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
1805:       fieldsplit preconditioner. If not set the matrix block size is used.

1807:     Logically Collective on PC

1809:     Input Parameters:
1810: +   pc  - the preconditioner context
1811: -   bs - the block size

1813:     Level: intermediate

1815: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()

1817: @*/
1818: PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1819: {

1825:   PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));
1826:   return(0);
1827: }

1831: /*@C
1832:    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits

1834:    Collective on KSP

1836:    Input Parameter:
1837: .  pc - the preconditioner context

1839:    Output Parameters:
1840: +  n - the number of splits
1841: -  pc - the array of KSP contexts

1843:    Note:
1844:    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1845:    (not the KSP just the array that contains them).

1847:    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().

1849:    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1850:       You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the
1851:       KSP array must be.


1854:    Level: advanced

1856: .seealso: PCFIELDSPLIT
1857: @*/
1858: PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1859: {

1865:   PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
1866:   return(0);
1867: }

1871: /*@
1872:     PCFieldSplitSetSchurPre -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1873:       A11 matrix. Otherwise no preconditioner is used.

1875:     Collective on PC

1877:     Input Parameters:
1878: +   pc      - the preconditioner context
1879: .   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_SCHUR_PRE_USER 
1880:               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
1881: -   userpre - matrix to use for preconditioning, or NULL

1883:     Options Database:
1884: .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments

1886:     Notes:
1887: $    If ptype is
1888: $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
1889: $             matrix associated with the Schur complement (i.e. A11), not he Schur complement matrix
1890: $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
1891: $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
1892: $             preconditioner
1893: $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
1894: $             to this function).
1895: $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
1896: $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
1897: $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
1898: $        full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PFIELDSPLIT (this is expensive)
1899: $             useful mostly as a test that the Schur complement approach can work for your problem

1901:      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1902:     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1903:     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.

1905:     Level: intermediate

1907: .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
1908:           MatSchurComplementSetAinvType(), PCLSC

1910: @*/
1911: PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1912: {

1917:   PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));
1918:   return(0);
1919: }
1920: PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */

1924: /*@
1925:     PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
1926:     preconditioned.  See PCFieldSplitSetSchurPre() for details.

1928:     Logically Collective on PC

1930:     Input Parameters:
1931: .   pc      - the preconditioner context

1933:     Output Parameters:
1934: +   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11, PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
1935: -   userpre - matrix to use for preconditioning (with PC_FIELDSPLIT_PRE_USER), or NULL

1937:     Level: intermediate

1939: .seealso: PCFieldSplitSetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC

1941: @*/
1942: PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
1943: {

1948:   PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));
1949:   return(0);
1950: }

1954: /*@
1955:     PCFieldSplitSchurGetS -  extract the MatSchurComplement object used by this PC in case it needs to be configured separately

1957:     Not collective

1959:     Input Parameter:
1960: .   pc      - the preconditioner context

1962:     Output Parameter:
1963: .   S       - the Schur complement matrix

1965:     Notes:
1966:     This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().

1968:     Level: advanced

1970: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurRestoreS()

1972: @*/
1973: PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
1974: {
1976:   const char*    t;
1977:   PetscBool      isfs;
1978:   PC_FieldSplit  *jac;

1982:   PetscObjectGetType((PetscObject)pc,&t);
1983:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
1984:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1985:   jac = (PC_FieldSplit*)pc->data;
1986:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1987:   if (S) *S = jac->schur;
1988:   return(0);
1989: }

1993: /*@
1994:     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC

1996:     Not collective

1998:     Input Parameters:
1999: +   pc      - the preconditioner context
2000: .   S       - the Schur complement matrix

2002:     Level: advanced

2004: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurGetS()

2006: @*/
2007: PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
2008: {
2010:   const char*    t;
2011:   PetscBool      isfs;
2012:   PC_FieldSplit  *jac;

2016:   PetscObjectGetType((PetscObject)pc,&t);
2017:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
2018:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2019:   jac = (PC_FieldSplit*)pc->data;
2020:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
2021:   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
2022:   return(0);
2023: }


2028: static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2029: {
2030:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2034:   jac->schurpre = ptype;
2035:   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
2036:     MatDestroy(&jac->schur_user);
2037:     jac->schur_user = pre;
2038:     PetscObjectReference((PetscObject)jac->schur_user);
2039:   }
2040:   return(0);
2041: }

2045: static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
2046: {
2047:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2050:   *ptype = jac->schurpre;
2051:   *pre   = jac->schur_user;
2052:   return(0);
2053: }

2057: /*@
2058:     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain

2060:     Collective on PC

2062:     Input Parameters:
2063: +   pc  - the preconditioner context
2064: -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default

2066:     Options Database:
2067: .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full


2070:     Level: intermediate

2072:     Notes:
2073:     The FULL factorization is

2075: $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
2076: $   (C   D)    (C*Ainv  1) (0   S) (0     1  )

2078:     where S = D - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D,
2079:     and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES).

2081:     If applied exactly, FULL factorization is a direct solver. The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial
2082:     of degree 2, so KSPGMRES converges in 2 iterations. If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
2083:     application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice. With DIAG,
2084:     the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.

2086:     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES. Note that a flexible method like KSPFGMRES
2087:     or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used inside a split).

2089:     References:
2090: +   1. - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
2091: -   2. - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).

2093: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
2094: @*/
2095: PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2096: {

2101:   PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));
2102:   return(0);
2103: }

2107: static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2108: {
2109:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2112:   jac->schurfactorization = ftype;
2113:   return(0);
2114: }

2118: /*@C
2119:    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement

2121:    Collective on KSP

2123:    Input Parameter:
2124: .  pc - the preconditioner context

2126:    Output Parameters:
2127: +  A00 - the (0,0) block
2128: .  A01 - the (0,1) block
2129: .  A10 - the (1,0) block
2130: -  A11 - the (1,1) block

2132:    Level: advanced

2134: .seealso: PCFIELDSPLIT
2135: @*/
2136: PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
2137: {
2138:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2142:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2143:   if (A00) *A00 = jac->pmat[0];
2144:   if (A01) *A01 = jac->B;
2145:   if (A10) *A10 = jac->C;
2146:   if (A11) *A11 = jac->pmat[1];
2147:   return(0);
2148: }

2152: static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
2153: {
2154:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2158:   jac->type = type;
2159:   if (type == PC_COMPOSITE_SCHUR) {
2160:     pc->ops->apply = PCApply_FieldSplit_Schur;
2161:     pc->ops->view  = PCView_FieldSplit_Schur;

2163:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);
2164:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);
2165:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);
2166:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);

2168:   } else {
2169:     pc->ops->apply = PCApply_FieldSplit;
2170:     pc->ops->view  = PCView_FieldSplit;

2172:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2173:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);
2174:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);
2175:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);
2176:   }
2177:   return(0);
2178: }

2182: static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
2183: {
2184:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2187:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2188:   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
2189:   jac->bs = bs;
2190:   return(0);
2191: }

2195: /*@
2196:    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.

2198:    Collective on PC

2200:    Input Parameter:
2201: .  pc - the preconditioner context
2202: .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

2204:    Options Database Key:
2205: .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type

2207:    Level: Intermediate

2209: .keywords: PC, set, type, composite preconditioner, additive, multiplicative

2211: .seealso: PCCompositeSetType()

2213: @*/
2214: PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
2215: {

2220:   PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));
2221:   return(0);
2222: }

2226: /*@
2227:   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.

2229:   Not collective

2231:   Input Parameter:
2232: . pc - the preconditioner context

2234:   Output Parameter:
2235: . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

2237:   Level: Intermediate

2239: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2240: .seealso: PCCompositeSetType()
2241: @*/
2242: PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2243: {
2244:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2249:   *type = jac->type;
2250:   return(0);
2251: }

2255: /*@
2256:    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.

2258:    Logically Collective

2260:    Input Parameters:
2261: +  pc   - the preconditioner context
2262: -  flg  - boolean indicating whether to use field splits defined by the DM

2264:    Options Database Key:
2265: .  -pc_fieldsplit_dm_splits

2267:    Level: Intermediate

2269: .keywords: PC, DM, composite preconditioner, additive, multiplicative

2271: .seealso: PCFieldSplitGetDMSplits()

2273: @*/
2274: PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
2275: {
2276:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2277:   PetscBool      isfs;

2283:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2284:   if (isfs) {
2285:     jac->dm_splits = flg;
2286:   }
2287:   return(0);
2288: }


2293: /*@
2294:    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.

2296:    Logically Collective

2298:    Input Parameter:
2299: .  pc   - the preconditioner context

2301:    Output Parameter:
2302: .  flg  - boolean indicating whether to use field splits defined by the DM

2304:    Level: Intermediate

2306: .keywords: PC, DM, composite preconditioner, additive, multiplicative

2308: .seealso: PCFieldSplitSetDMSplits()

2310: @*/
2311: PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
2312: {
2313:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2314:   PetscBool      isfs;

2320:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2321:   if (isfs) {
2322:     if(flg) *flg = jac->dm_splits;
2323:   }
2324:   return(0);
2325: }

2327: /* -------------------------------------------------------------------------------------*/
2328: /*MC
2329:    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2330:                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.

2332:      To set options on the solvers for each block append -fieldsplit_ to all the PC
2333:         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1

2335:      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
2336:          and set the options directly on the resulting KSP object

2338:    Level: intermediate

2340:    Options Database Keys:
2341: +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
2342: .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
2343:                               been supplied explicitly by -pc_fieldsplit_%d_fields
2344: .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
2345: .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
2346: .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see PCFieldSplitSetSchurPre()
2347: .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver

2349: -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
2350:      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields

2352:    Notes:
2353:     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
2354:      to define a field by an arbitrary collection of entries.

2356:       If no fields are set the default is used. The fields are defined by entries strided by bs,
2357:       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
2358:       if this is not called the block size defaults to the blocksize of the second matrix passed
2359:       to KSPSetOperators()/PCSetOperators().

2361: $     For the Schur complement preconditioner if J = ( A00 A01 )
2362: $                                                    ( A10 A11 )
2363: $     the preconditioner using full factorization is
2364: $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2365: $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
2366:      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
2367: $              S = A11 - A10 ksp(A00) A01
2368:      which is usually dense and not stored explicitly.  The action of ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
2369:      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2370:      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2371:      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.

2373:      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
2374:      diag gives
2375: $              ( inv(A00)     0   )
2376: $              (   0      -ksp(S) )
2377:      note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of
2378: $              (  A00   0 )
2379: $              (  A10   S )
2380:      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2381: $              ( A00 A01 )
2382: $              (  0   S  )
2383:      where again the inverses of A00 and S are applied using KSPs.

2385:      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
2386:      is used automatically for a second block.

2388:      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
2389:      Generally it should be used with the AIJ format.

2391:      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
2392:      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
2393:      inside a smoother resulting in "Distributive Smoothers".

2395:    Concepts: physics based preconditioners, block preconditioners

2397:    There is a nice discussion of block preconditioners in

2399: [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
2400:        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
2401:        http://chess.cs.umd.edu/~elman/papers/tax.pdf

2403:    The Constrained Pressure Preconditioner (CPR) does not appear to be currently implementable directly with PCFIELDSPLIT. CPR solves first the Schur complemented pressure equation, updates the
2404:    residual on all variables and then applies a simple ILU like preconditioner on all the variables. So it is very much like the full Schur complement with selfp representing the Schur complement but instead
2405:    of backsolving for the saturations in the last step it solves a full coupled (ILU) system for updates to all the variables.

2407: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
2408:            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
2409:            MatSchurComplementSetAinvType()
2410: M*/

2414: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
2415: {
2417:   PC_FieldSplit  *jac;

2420:   PetscNewLog(pc,&jac);

2422:   jac->bs                 = -1;
2423:   jac->nsplits            = 0;
2424:   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2425:   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2426:   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2427:   jac->dm_splits          = PETSC_TRUE;

2429:   pc->data = (void*)jac;

2431:   pc->ops->apply           = PCApply_FieldSplit;
2432:   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
2433:   pc->ops->setup           = PCSetUp_FieldSplit;
2434:   pc->ops->reset           = PCReset_FieldSplit;
2435:   pc->ops->destroy         = PCDestroy_FieldSplit;
2436:   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
2437:   pc->ops->view            = PCView_FieldSplit;
2438:   pc->ops->applyrichardson = 0;

2440:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2441:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);
2442:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);
2443:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);
2444:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);
2445:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);
2446:   return(0);
2447: }