Actual source code: fieldsplit.c

petsc-master 2016-05-28
Report Typos and Errors
  2: #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
  3: #include <petsc/private/kspimpl.h>
  4: #include <petscdm.h>


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

 10: 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;

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

496:     jac->issetup = PETSC_TRUE;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

713:       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
714:       MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);
715:       MatSetType(jac->schur,MATSCHURCOMPLEMENT);
716:       MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
717:       PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
718:       /* 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? */
719:       MatSetOptionsPrefix(jac->schur,schurmatprefix);
720:       MatSetFromOptions(jac->schur);
721:       MatGetNullSpace(jac->mat[1], &sp);
722:       if (sp) {
723:         MatSetNullSpace(jac->schur, sp);
724:       }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1258:   PetscOptionsHead(PetscOptionsObject,"FieldSplit options");
1259:   PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);
1260:   PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);
1261:   if (flg) {
1262:     PCFieldSplitSetBlockSize(pc,bs);
1263:   }
1264:   jac->diag_use_amat = pc->useAmat;
1265:   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);
1266:   jac->offdiag_use_amat = pc->useAmat;
1267:   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);
1268:   /* FIXME: No programmatic equivalent to the following. */
1269:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);
1270:   if (stokes) {
1271:     PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
1272:     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1273:   }

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

1296: /*------------------------------------------------------------------------------------*/

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

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

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

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

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

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

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

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

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

1381:   PetscMalloc1(jac->nsplits,subksp);
1382:   while (ilink) {
1383:     (*subksp)[cnt++] = ilink->ksp;
1384:     ilink            = ilink->next;
1385:   }
1386:   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);
1387:   if (n) *n = jac->nsplits;
1388:   return(0);
1389: }

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

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

1400:     Level: advanced

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

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


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

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

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

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

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

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

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

1545:     Logically Collective on PC

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

1553:     Level: intermediate

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

1557:      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1558:      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
1559:      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1560:      where the numbered entries indicate what is in the field.

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

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

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

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

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

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

1590:     Logically Collective on PC

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

1596:     Options Database:
1597: .     -pc_fieldsplit_diag_use_amat

1599:     Level: intermediate

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

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

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

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

1623:     Logically Collective on PC

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

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


1632:     Level: intermediate

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

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

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

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

1657:     Logically Collective on PC

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

1663:     Options Database:
1664: .     -pc_fieldsplit_off_diag_use_amat

1666:     Level: intermediate

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

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

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

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

1690:     Logically Collective on PC

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

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


1699:     Level: intermediate

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

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

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



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

1726:     Logically Collective on PC

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


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

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

1740:     Level: intermediate

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

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

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

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

1762:     Logically Collective on PC

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

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

1771:     Level: intermediate

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

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

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

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

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

1808:     Logically Collective on PC

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

1814:     Level: intermediate

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

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

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

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

1835:    Collective on KSP

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

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

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

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

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


1855:    Level: advanced

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

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

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

1876:     Collective on PC

1878:     Input Parameters:
1879: +   pc      - the preconditioner context
1880: .   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
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

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

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

1904:     Level: intermediate

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

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

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

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

1927:     Logically Collective on PC

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

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

1936:     Level: intermediate

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

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

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

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

1956:     Not collective

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

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

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

1967:     Level: advanced

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

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

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

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

1995:     Not collective

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

2001:     Level: advanced

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

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

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


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

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

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

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

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

2059:     Collective on PC

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

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


2069:     Level: intermediate

2071:     Notes:
2072:     The FULL factorization is

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

2077:     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,
2078:     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).

2080:     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
2081:     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
2082:     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,
2083:     the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.

2085:     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
2086:     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).

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

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

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

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

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

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

2120:    Collective on KSP

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

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

2131:    Level: advanced

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

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

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

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

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

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

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

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

2186:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2187:   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);
2188:   jac->bs = bs;
2189:   return(0);
2190: }

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

2197:    Collective on PC

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

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

2206:    Level: Intermediate

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

2210: .seealso: PCCompositeSetType()

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

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

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

2228:   Not collective

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

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

2236:   Level: Intermediate

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

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

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

2257:    Logically Collective

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

2263:    Options Database Key:
2264: .  -pc_fieldsplit_dm_splits

2266:    Level: Intermediate

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

2270: .seealso: PCFieldSplitGetDMSplits()

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

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


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

2295:    Logically Collective

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

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

2303:    Level: Intermediate

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

2307: .seealso: PCFieldSplitSetDMSplits()

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

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

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

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

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

2337:    Level: intermediate

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

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

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

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

2360: $     For the Schur complement preconditioner if J = ( A00 A01 )
2361: $                                                    ( A10 A11 )
2362: $     the preconditioner using full factorization is
2363: $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2364: $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
2365:      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
2366: $              S = A11 - A10 ksp(A00) A01
2367:      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
2368:      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2369:      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2370:      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.

2372:      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
2373:      diag gives
2374: $              ( inv(A00)     0   )
2375: $              (   0      -ksp(S) )
2376:      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
2377: $              (  A00   0 )
2378: $              (  A10   S )
2379:      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2380: $              ( A00 A01 )
2381: $              (  0   S  )
2382:      where again the inverses of A00 and S are applied using KSPs.

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

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

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

2394:    Concepts: physics based preconditioners, block preconditioners

2396:    There is a nice discussion of block preconditioners in

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

2402:    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
2403:    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
2404:    of backsolving for the saturations in the last step it solves a full coupled (ILU) system for updates to all the variables.

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

2413: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
2414: {
2416:   PC_FieldSplit  *jac;

2419:   PetscNewLog(pc,&jac);

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

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

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

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