Actual source code: fieldsplit.c

petsc-3.10.2 2018-10-09
Report Typos and Errors


  3:  #include <petsc/private/pcimpl.h>
  4: #include <petsc/private/kspimpl.h>    /*  This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
  5:  #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;
 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:   PetscScalar               schurscale;            /* Scaling factor for the Schur complement solution with DIAG factorization */

 49:   PC_FieldSplitLink         head;
 50:   PetscBool                 isrestrict;             /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
 51:   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
 52:   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
 53:   PetscBool                 diag_use_amat;          /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
 54:   PetscBool                 offdiag_use_amat;       /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
 55:   PetscBool                 detect;                 /* Whether to form 2-way split by finding zero diagonal entries */
 56: } PC_FieldSplit;

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

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


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

 91:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 92:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
 93:   if (iascii) {
 94:     if (jac->bs > 0) {
 95:       PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);
 96:     } else {
 97:       PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);
 98:     }
 99:     if (pc->useAmat) {
100:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");
101:     }
102:     if (jac->diag_use_amat) {
103:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for diagonal blocks\n");
104:     }
105:     if (jac->offdiag_use_amat) {
106:       PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for off-diagonal blocks\n");
107:     }
108:     PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");
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:   }

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

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

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

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

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

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

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

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: }

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

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

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

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

401:       if (jac->detect) {
402:         IS       zerodiags,rest;
403:         PetscInt nmin,nmax;

405:         MatGetOwnershipRange(pc->mat,&nmin,&nmax);
406:         MatFindZeroDiagonals(pc->mat,&zerodiags);
407:         ISComplement(zerodiags,nmin,nmax,&rest);
408:         PCFieldSplitSetIS(pc,"0",rest);
409:         PCFieldSplitSetIS(pc,"1",zerodiags);
410:         ISDestroy(&zerodiags);
411:         ISDestroy(&rest);
412:       } else if (coupling) {
413:         IS       coupling,rest;
414:         PetscInt nmin,nmax;

416:         MatGetOwnershipRange(pc->mat,&nmin,&nmax);
417:         MatFindOffBlockDiagonalEntries(pc->mat,&coupling);
418:         ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);
419:         ISSetIdentity(rest);
420:         PCFieldSplitSetIS(pc,"0",rest);
421:         PCFieldSplitSetIS(pc,"1",coupling);
422:         ISDestroy(&coupling);
423:         ISDestroy(&rest);
424:       } else {
425:         PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);
426:         if (!fieldsplit_default) {
427:           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
428:            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
429:           PCFieldSplitSetRuntimeSplits_Private(pc);
430:           if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
431:         }
432:         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
433:           PetscInfo(pc,"Using default splitting of fields\n");
434:           for (i=0; i<jac->bs; i++) {
435:             char splitname[8];
436:             PetscSNPrintf(splitname,sizeof(splitname),"%D",i);
437:             PCFieldSplitSetFields(pc,splitname,1,&i,&i);
438:           }
439:           jac->defaultsplit = PETSC_TRUE;
440:         }
441:       }
442:     }
443:   } else if (jac->nsplits == 1) {
444:     if (ilink->is) {
445:       IS       is2;
446:       PetscInt nmin,nmax;

448:       MatGetOwnershipRange(pc->mat,&nmin,&nmax);
449:       ISComplement(ilink->is,nmin,nmax,&is2);
450:       PCFieldSplitSetIS(pc,"1",is2);
451:       ISDestroy(&is2);
452:     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
453:   }

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

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

461: static PetscErrorCode PCSetUp_FieldSplit(PC pc)
462: {
463:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
464:   PetscErrorCode    ierr;
465:   PC_FieldSplitLink ilink;
466:   PetscInt          i,nsplit;
467:   PetscBool         sorted, sorted_col;

470:   PCFieldSplitSetDefaults(pc);
471:   nsplit = jac->nsplits;
472:   ilink  = jac->head;

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

478:     jac->issetup = PETSC_TRUE;

480:     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
481:     if (jac->defaultsplit || !ilink->is) {
482:       if (jac->bs <= 0) jac->bs = nsplit;
483:     }
484:     bs     = jac->bs;
485:     MatGetOwnershipRange(pc->pmat,&rstart,&rend);
486:     nslots = (rend - rstart)/bs;
487:     for (i=0; i<nsplit; i++) {
488:       if (jac->defaultsplit) {
489:         ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);
490:         ISDuplicate(ilink->is,&ilink->is_col);
491:       } else if (!ilink->is) {
492:         if (ilink->nfields > 1) {
493:           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
494:           PetscMalloc1(ilink->nfields*nslots,&ii);
495:           PetscMalloc1(ilink->nfields*nslots,&jj);
496:           for (j=0; j<nslots; j++) {
497:             for (k=0; k<nfields; k++) {
498:               ii[nfields*j + k] = rstart + bs*j + fields[k];
499:               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
500:             }
501:           }
502:           ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);
503:           ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);
504:           ISSetBlockSize(ilink->is, nfields);
505:           ISSetBlockSize(ilink->is_col, nfields);
506:         } else {
507:           ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);
508:           ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);
509:        }
510:       }
511:       ISSorted(ilink->is,&sorted);
512:       if (ilink->is_col) { ISSorted(ilink->is_col,&sorted_col); }
513:       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
514:       ilink = ilink->next;
515:     }
516:   }

518:   ilink = jac->head;
519:   if (!jac->pmat) {
520:     Vec xtmp;

522:     MatCreateVecs(pc->pmat,&xtmp,NULL);
523:     PetscMalloc1(nsplit,&jac->pmat);
524:     PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);
525:     for (i=0; i<nsplit; i++) {
526:       MatNullSpace sp;

528:       /* Check for preconditioning matrix attached to IS */
529:       PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);
530:       if (jac->pmat[i]) {
531:         PetscObjectReference((PetscObject) jac->pmat[i]);
532:         if (jac->type == PC_COMPOSITE_SCHUR) {
533:           jac->schur_user = jac->pmat[i];

535:           PetscObjectReference((PetscObject) jac->schur_user);
536:         }
537:       } else {
538:         const char *prefix;
539:         MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);
540:         KSPGetOptionsPrefix(ilink->ksp,&prefix);
541:         MatSetOptionsPrefix(jac->pmat[i],prefix);
542:         MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");
543:       }
544:       /* create work vectors for each split */
545:       MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);
546:       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
547:       /* compute scatter contexts needed by multiplicative versions and non-default splits */
548:       VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);
549:       PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);
550:       if (sp) {
551:         MatSetNearNullSpace(jac->pmat[i], sp);
552:       }
553:       ilink = ilink->next;
554:     }
555:     VecDestroy(&xtmp);
556:   } else {
557:     MatReuse scall;
558:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
559:       for (i=0; i<nsplit; i++) {
560:         MatDestroy(&jac->pmat[i]);
561:       }
562:       scall = MAT_INITIAL_MATRIX;
563:     } else scall = MAT_REUSE_MATRIX;

565:     for (i=0; i<nsplit; i++) {
566:       Mat pmat;

568:       /* Check for preconditioning matrix attached to IS */
569:       PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);
570:       if (!pmat) {
571:         MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,scall,&jac->pmat[i]);
572:       }
573:       ilink = ilink->next;
574:     }
575:   }
576:   if (jac->diag_use_amat) {
577:     ilink = jac->head;
578:     if (!jac->mat) {
579:       PetscMalloc1(nsplit,&jac->mat);
580:       for (i=0; i<nsplit; i++) {
581:         MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);
582:         ilink = ilink->next;
583:       }
584:     } else {
585:       MatReuse scall;
586:       if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
587:         for (i=0; i<nsplit; i++) {
588:           MatDestroy(&jac->mat[i]);
589:         }
590:         scall = MAT_INITIAL_MATRIX;
591:       } else scall = MAT_REUSE_MATRIX;

593:       for (i=0; i<nsplit; i++) {
594:         if (jac->mat[i]) {MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,scall,&jac->mat[i]);}
595:         ilink = ilink->next;
596:       }
597:     }
598:   } else {
599:     jac->mat = jac->pmat;
600:   }

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

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

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

636:         if (jac->offdiag_use_amat) {
637:           MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,scall,&jac->Afield[1]);
638:         } else {
639:           MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,scall,&jac->Afield[1]);
640:         }
641:       }
642:     } else {
643:       if (!jac->Afield) {
644:         PetscMalloc1(nsplit,&jac->Afield);
645:         for (i=0; i<nsplit; i++) {
646:           if (jac->offdiag_use_amat) {
647:             MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);
648:           } else {
649:             MatCreateSubMatrix(pc->pmat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);
650:           }
651:           ilink = ilink->next;
652:         }
653:       } else {
654:         MatReuse scall;
655:         if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
656:           for (i=0; i<nsplit; i++) {
657:             MatDestroy(&jac->Afield[i]);
658:           }
659:           scall = MAT_INITIAL_MATRIX;
660:         } else scall = MAT_REUSE_MATRIX;

662:         for (i=0; i<nsplit; i++) {
663:           if (jac->offdiag_use_amat) {
664:             MatCreateSubMatrix(pc->mat,ilink->is,NULL,scall,&jac->Afield[i]);
665:           } else {
666:             MatCreateSubMatrix(pc->pmat,ilink->is,NULL,scall,&jac->Afield[i]);
667:           }
668:           ilink = ilink->next;
669:         }
670:       }
671:     }
672:   }

674:   if (jac->type == PC_COMPOSITE_SCHUR) {
675:     IS          ccis;
676:     PetscBool   isspd;
677:     PetscInt    rstart,rend;
678:     char        lscname[256];
679:     PetscObject LSC_L;

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

683:     /* If pc->mat is SPD, don't scale by -1 the Schur complement */
684:     if (jac->schurscale == (PetscScalar)-1.0) {
685:       MatGetOption(pc->pmat,MAT_SPD,&isspd);
686:       jac->schurscale = (isspd == PETSC_TRUE) ? 1.0 : -1.0;
687:     }

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

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

696:       MatSchurComplementGetKSP(jac->schur, &kspInner);
697:       ilink = jac->head;
698:       ISComplement(ilink->is_col,rstart,rend,&ccis);
699:       if (jac->offdiag_use_amat) {
700:         MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);
701:       } else {
702:         MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);
703:       }
704:       ISDestroy(&ccis);
705:       ilink = ilink->next;
706:       ISComplement(ilink->is_col,rstart,rend,&ccis);
707:       if (jac->offdiag_use_amat) {
708:         MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);
709:       } else {
710:         MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);
711:       }
712:       ISDestroy(&ccis);
713:       MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
714:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
715:         MatDestroy(&jac->schurp);
716:         MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);
717:       }
718:       if (kspA != kspInner) {
719:         KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);
720:       }
721:       if (kspUpper != kspA) {
722:         KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);
723:       }
724:       KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));
725:     } else {
726:       const char   *Dprefix;
727:       char         schurprefix[256], schurmatprefix[256];
728:       char         schurtestoption[256];
729:       MatNullSpace sp;
730:       PetscBool    flg;
731:       KSP          kspt;

733:       /* extract the A01 and A10 matrices */
734:       ilink = jac->head;
735:       ISComplement(ilink->is_col,rstart,rend,&ccis);
736:       if (jac->offdiag_use_amat) {
737:         MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
738:       } else {
739:         MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
740:       }
741:       ISDestroy(&ccis);
742:       ilink = ilink->next;
743:       ISComplement(ilink->is_col,rstart,rend,&ccis);
744:       if (jac->offdiag_use_amat) {
745:         MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
746:       } else {
747:         MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
748:       }
749:       ISDestroy(&ccis);

751:       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
752:       MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);
753:       MatSetType(jac->schur,MATSCHURCOMPLEMENT);
754:       MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
755:       PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
756:       MatSetOptionsPrefix(jac->schur,schurmatprefix);
757:       MatSchurComplementGetKSP(jac->schur,&kspt);
758:       KSPSetOptionsPrefix(kspt,schurmatprefix);

760:       /* Note: this is not true in general */
761:       MatGetNullSpace(jac->mat[1], &sp);
762:       if (sp) {
763:         MatSetNullSpace(jac->schur, sp);
764:       }

766:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);
767:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
768:       if (flg) {
769:         DM  dmInner;
770:         KSP kspInner;
771:         PC  pcInner;

773:         MatSchurComplementGetKSP(jac->schur, &kspInner);
774:         KSPReset(kspInner);
775:         KSPSetOperators(kspInner,jac->mat[0],jac->pmat[0]);
776:         PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
777:         /* Indent this deeper to emphasize the "inner" nature of this solver. */
778:         PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);
779:         PetscObjectIncrementTabLevel((PetscObject)kspInner->pc, (PetscObject) pc, 2);
780:         KSPSetOptionsPrefix(kspInner, schurprefix);

782:         /* Set DM for new solver */
783:         KSPGetDM(jac->head->ksp, &dmInner);
784:         KSPSetDM(kspInner, dmInner);
785:         KSPSetDMActive(kspInner, PETSC_FALSE);

787:         /* Defaults to PCKSP as preconditioner */
788:         KSPGetPC(kspInner, &pcInner);
789:         PCSetType(pcInner, PCKSP);
790:         PCKSPSetKSP(pcInner, jac->head->ksp);
791:       } else {
792:          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
793:           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
794:           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
795:           * S = D - C A_inner^{-1} B, we expect S to be defined using an accurate definition of A_inner^{-1}, so we make
796:           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
797:           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
798:         KSPSetType(jac->head->ksp,KSPGMRES);
799:         MatSchurComplementSetKSP(jac->schur,jac->head->ksp);
800:       }
801:       KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);
802:       KSPSetFromOptions(jac->head->ksp);
803:       MatSetFromOptions(jac->schur);

805:       PetscObjectTypeCompare((PetscObject)jac->schur, MATSCHURCOMPLEMENT, &flg);
806:       if (flg) { /* Need to do this otherwise PCSetUp_KSP will overwrite the amat of jac->head->ksp */
807:         KSP kspInner;
808:         PC  pcInner;

810:         MatSchurComplementGetKSP(jac->schur, &kspInner);
811:         KSPGetPC(kspInner, &pcInner);
812:         PetscObjectTypeCompare((PetscObject)pcInner, PCKSP, &flg);
813:         if (flg) {
814:           KSP ksp;

816:           PCKSPGetKSP(pcInner, &ksp);
817:           if (ksp == jac->head->ksp) {
818:             PCSetUseAmat(pcInner, PETSC_TRUE);
819:           }
820:         }
821:       }
822:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);
823:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
824:       if (flg) {
825:         DM dmInner;

827:         PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
828:         KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);
829:         KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);
830:         KSPSetOptionsPrefix(jac->kspupper, schurprefix);
831:         PetscObjectIncrementTabLevel((PetscObject)jac->kspupper, (PetscObject) pc, 1);
832:         PetscObjectIncrementTabLevel((PetscObject)jac->kspupper->pc, (PetscObject) pc, 1);
833:         KSPGetDM(jac->head->ksp, &dmInner);
834:         KSPSetDM(jac->kspupper, dmInner);
835:         KSPSetDMActive(jac->kspupper, PETSC_FALSE);
836:         KSPSetFromOptions(jac->kspupper);
837:         KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);
838:         VecDuplicate(jac->head->x, &jac->head->z);
839:       } else {
840:         jac->kspupper = jac->head->ksp;
841:         PetscObjectReference((PetscObject) jac->head->ksp);
842:       }

844:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
845:         MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);
846:       }
847:       KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);
848:       KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);
849:       PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);
850:       PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);
851:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
852:         PC pcschur;
853:         KSPGetPC(jac->kspschur,&pcschur);
854:         PCSetType(pcschur,PCNONE);
855:         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
856:       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
857:         MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);
858:       }
859:       KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));
860:       KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);
861:       KSPSetOptionsPrefix(jac->kspschur,         Dprefix);
862:       /* propagate DM */
863:       {
864:         DM sdm;
865:         KSPGetDM(jac->head->next->ksp, &sdm);
866:         if (sdm) {
867:           KSPSetDM(jac->kspschur, sdm);
868:           KSPSetDMActive(jac->kspschur, PETSC_FALSE);
869:         }
870:       }
871:       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
872:       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
873:       KSPSetFromOptions(jac->kspschur);
874:     }
875:     MatAssemblyBegin(jac->schur,MAT_FINAL_ASSEMBLY);
876:     MatAssemblyEnd(jac->schur,MAT_FINAL_ASSEMBLY);

878:     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
879:     PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);
880:     PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);
881:     if (!LSC_L) {PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);}
882:     if (LSC_L) {PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);}
883:     PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);
884:     PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);
885:     if (!LSC_L) {PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);}
886:     if (LSC_L) {PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);}
887:   } else {
888:     /* set up the individual splits' PCs */
889:     i     = 0;
890:     ilink = jac->head;
891:     while (ilink) {
892:       KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);
893:       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
894:       if (!jac->suboptionsset) {KSPSetFromOptions(ilink->ksp);}
895:       i++;
896:       ilink = ilink->next;
897:     }
898:   }

900:   jac->suboptionsset = PETSC_TRUE;
901:   return(0);
902: }

904: #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
905:   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
906:    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
907:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
908:    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
909:    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
910:    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
911:    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))

913: static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
914: {
915:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
916:   PetscErrorCode    ierr;
917:   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
918:   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;

921:   switch (jac->schurfactorization) {
922:   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
923:     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
924:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
925:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
926:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
927:     PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
928:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
929:     PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
930:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
931:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
932:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
933:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
934:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
935:     VecScale(ilinkD->y,jac->schurscale);
936:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
937:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
938:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
939:     break;
940:   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
941:     /* [A00 0; A10 S], suitable for left preconditioning */
942:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
943:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
944:     PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
945:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
946:     PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
947:     MatMult(jac->C,ilinkA->y,ilinkD->x);
948:     VecScale(ilinkD->x,-1.);
949:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
950:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
951:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
952:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
953:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
954:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
955:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
956:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
957:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
958:     break;
959:   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
960:     /* [A00 A01; 0 S], suitable for right preconditioning */
961:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
962:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
963:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
964:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
965:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);    MatMult(jac->B,ilinkD->y,ilinkA->x);
966:     VecScale(ilinkA->x,-1.);
967:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
968:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
969:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
970:     PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
971:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
972:     PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
973:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
974:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
975:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
976:     break;
977:   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
978:     /* [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 */
979:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
980:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
981:     PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);
982:     KSPSolve(kspLower,ilinkA->x,ilinkA->y);
983:     PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);
984:     MatMult(jac->C,ilinkA->y,ilinkD->x);
985:     VecScale(ilinkD->x,-1.0);
986:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
987:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);

989:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
990:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
991:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
992:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
993:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);

995:     if (kspUpper == kspA) {
996:       MatMult(jac->B,ilinkD->y,ilinkA->y);
997:       VecAXPY(ilinkA->x,-1.0,ilinkA->y);
998:       PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
999:       KSPSolve(kspA,ilinkA->x,ilinkA->y);
1000:       PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1001:     } else {
1002:       PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1003:       KSPSolve(kspA,ilinkA->x,ilinkA->y);
1004:       PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1005:       MatMult(jac->B,ilinkD->y,ilinkA->x);
1006:       PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);
1007:       KSPSolve(kspUpper,ilinkA->x,ilinkA->z);
1008:       PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);
1009:       VecAXPY(ilinkA->y,-1.0,ilinkA->z);
1010:     }
1011:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1012:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1013:   }
1014:   return(0);
1015: }

1017: static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
1018: {
1019:   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1020:   PetscErrorCode     ierr;
1021:   PC_FieldSplitLink  ilink = jac->head;
1022:   PetscInt           cnt,bs;
1023:   KSPConvergedReason reason;

1026:   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1027:     if (jac->defaultsplit) {
1028:       VecGetBlockSize(x,&bs);
1029:       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);
1030:       VecGetBlockSize(y,&bs);
1031:       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);
1032:       VecStrideGatherAll(x,jac->x,INSERT_VALUES);
1033:       while (ilink) {
1034:         PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1035:         KSPSolve(ilink->ksp,ilink->x,ilink->y);
1036:         PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1037:         KSPGetConvergedReason(ilink->ksp,&reason);
1038:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1039:           pc->failedreason = PC_SUBPC_ERROR;
1040:         }
1041:         ilink = ilink->next;
1042:       }
1043:       VecStrideScatterAll(jac->y,y,INSERT_VALUES);
1044:     } else {
1045:       VecSet(y,0.0);
1046:       while (ilink) {
1047:         FieldSplitSplitSolveAdd(ilink,x,y);
1048:         KSPGetConvergedReason(ilink->ksp,&reason);
1049:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1050:           pc->failedreason = PC_SUBPC_ERROR;
1051:         }
1052:         ilink = ilink->next;
1053:       }
1054:     }
1055:   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
1056:     VecSet(y,0.0);
1057:     /* solve on first block for first block variables */
1058:     VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
1059:     VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
1060:     PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1061:     KSPSolve(ilink->ksp,ilink->x,ilink->y);
1062:     PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1063:     KSPGetConvergedReason(ilink->ksp,&reason);
1064:     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1065:       pc->failedreason = PC_SUBPC_ERROR;
1066:     }
1067:     VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1068:     VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);

1070:     /* compute the residual only onto second block variables using first block variables */
1071:     MatMult(jac->Afield[1],ilink->y,ilink->next->x);
1072:     ilink = ilink->next;
1073:     VecScale(ilink->x,-1.0);
1074:     VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1075:     VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);

1077:     /* solve on second block variables */
1078:     PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1079:     KSPSolve(ilink->ksp,ilink->x,ilink->y);
1080:     PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1081:     KSPGetConvergedReason(ilink->ksp,&reason);
1082:     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1083:       pc->failedreason = PC_SUBPC_ERROR;
1084:     }
1085:     VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1086:     VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1087:   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1088:     if (!jac->w1) {
1089:       VecDuplicate(x,&jac->w1);
1090:       VecDuplicate(x,&jac->w2);
1091:     }
1092:     VecSet(y,0.0);
1093:     FieldSplitSplitSolveAdd(ilink,x,y);
1094:     KSPGetConvergedReason(ilink->ksp,&reason);
1095:     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1096:       pc->failedreason = PC_SUBPC_ERROR;
1097:     }
1098:     cnt  = 1;
1099:     while (ilink->next) {
1100:       ilink = ilink->next;
1101:       /* compute the residual only over the part of the vector needed */
1102:       MatMult(jac->Afield[cnt++],y,ilink->x);
1103:       VecScale(ilink->x,-1.0);
1104:       VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1105:       VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1106:       PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1107:       KSPSolve(ilink->ksp,ilink->x,ilink->y);
1108:       PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1109:       KSPGetConvergedReason(ilink->ksp,&reason);
1110:       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1111:         pc->failedreason = PC_SUBPC_ERROR;
1112:       }
1113:       VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1114:       VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1115:     }
1116:     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1117:       cnt -= 2;
1118:       while (ilink->previous) {
1119:         ilink = ilink->previous;
1120:         /* compute the residual only over the part of the vector needed */
1121:         MatMult(jac->Afield[cnt--],y,ilink->x);
1122:         VecScale(ilink->x,-1.0);
1123:         VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1124:         VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1125:         PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1126:         KSPSolve(ilink->ksp,ilink->x,ilink->y);
1127:         PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1128:         KSPGetConvergedReason(ilink->ksp,&reason);
1129:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1130:           pc->failedreason = PC_SUBPC_ERROR;
1131:         }
1132:         VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1133:         VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1134:       }
1135:     }
1136:   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
1137:   return(0);
1138: }

1140: #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1141:   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1142:    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1143:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1144:    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) ||                  \
1145:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1146:    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1147:    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))

1149: static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1150: {
1151:   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1152:   PetscErrorCode     ierr;
1153:   PC_FieldSplitLink  ilink = jac->head;
1154:   PetscInt           bs;
1155:   KSPConvergedReason reason;

1158:   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1159:     if (jac->defaultsplit) {
1160:       VecGetBlockSize(x,&bs);
1161:       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);
1162:       VecGetBlockSize(y,&bs);
1163:       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);
1164:       VecStrideGatherAll(x,jac->x,INSERT_VALUES);
1165:       while (ilink) {
1166:         PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1167:         KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);
1168:         PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1169:         KSPGetConvergedReason(ilink->ksp,&reason);
1170:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1171:           pc->failedreason = PC_SUBPC_ERROR;
1172:         }
1173:         ilink = ilink->next;
1174:       }
1175:       VecStrideScatterAll(jac->y,y,INSERT_VALUES);
1176:     } else {
1177:       VecSet(y,0.0);
1178:       while (ilink) {
1179:         FieldSplitSplitSolveAddTranspose(ilink,x,y);
1180:         KSPGetConvergedReason(ilink->ksp,&reason);
1181:         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1182:           pc->failedreason = PC_SUBPC_ERROR;
1183:         }
1184:         ilink = ilink->next;
1185:       }
1186:     }
1187:   } else {
1188:     if (!jac->w1) {
1189:       VecDuplicate(x,&jac->w1);
1190:       VecDuplicate(x,&jac->w2);
1191:     }
1192:     VecSet(y,0.0);
1193:     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1194:       FieldSplitSplitSolveAddTranspose(ilink,x,y);
1195:       KSPGetConvergedReason(ilink->ksp,&reason);
1196:       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1197:         pc->failedreason = PC_SUBPC_ERROR;
1198:       }
1199:       while (ilink->next) {
1200:         ilink = ilink->next;
1201:         MatMultTranspose(pc->mat,y,jac->w1);
1202:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1203:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1204:       }
1205:       while (ilink->previous) {
1206:         ilink = ilink->previous;
1207:         MatMultTranspose(pc->mat,y,jac->w1);
1208:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1209:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1210:       }
1211:     } else {
1212:       while (ilink->next) {   /* get to last entry in linked list */
1213:         ilink = ilink->next;
1214:       }
1215:       FieldSplitSplitSolveAddTranspose(ilink,x,y);
1216:       KSPGetConvergedReason(ilink->ksp,&reason);
1217:       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1218:         pc->failedreason = PC_SUBPC_ERROR;
1219:       }
1220:       while (ilink->previous) {
1221:         ilink = ilink->previous;
1222:         MatMultTranspose(pc->mat,y,jac->w1);
1223:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1224:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1225:       }
1226:     }
1227:   }
1228:   return(0);
1229: }

1231: static PetscErrorCode PCReset_FieldSplit(PC pc)
1232: {
1233:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1234:   PetscErrorCode    ierr;
1235:   PC_FieldSplitLink ilink = jac->head,next;

1238:   while (ilink) {
1239:     KSPDestroy(&ilink->ksp);
1240:     VecDestroy(&ilink->x);
1241:     VecDestroy(&ilink->y);
1242:     VecDestroy(&ilink->z);
1243:     VecScatterDestroy(&ilink->sctx);
1244:     ISDestroy(&ilink->is);
1245:     ISDestroy(&ilink->is_col);
1246:     PetscFree(ilink->splitname);
1247:     PetscFree(ilink->fields);
1248:     PetscFree(ilink->fields_col);
1249:     next  = ilink->next;
1250:     PetscFree(ilink);
1251:     ilink = next;
1252:   }
1253:   jac->head = NULL;
1254:   PetscFree2(jac->x,jac->y);
1255:   if (jac->mat && jac->mat != jac->pmat) {
1256:     MatDestroyMatrices(jac->nsplits,&jac->mat);
1257:   } else if (jac->mat) {
1258:     jac->mat = NULL;
1259:   }
1260:   if (jac->pmat) {MatDestroyMatrices(jac->nsplits,&jac->pmat);}
1261:   if (jac->Afield) {MatDestroyMatrices(jac->nsplits,&jac->Afield);}
1262:   jac->nsplits = 0;
1263:   VecDestroy(&jac->w1);
1264:   VecDestroy(&jac->w2);
1265:   MatDestroy(&jac->schur);
1266:   MatDestroy(&jac->schurp);
1267:   MatDestroy(&jac->schur_user);
1268:   KSPDestroy(&jac->kspschur);
1269:   KSPDestroy(&jac->kspupper);
1270:   MatDestroy(&jac->B);
1271:   MatDestroy(&jac->C);
1272:   jac->isrestrict = PETSC_FALSE;
1273:   return(0);
1274: }

1276: static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1277: {
1278:   PetscErrorCode    ierr;

1281:   PCReset_FieldSplit(pc);
1282:   PetscFree(pc->data);
1283:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurGetSubKSP_C",NULL);
1284:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);
1285:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);
1286:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);
1287:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);
1288:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);
1289:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);
1290:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);
1291:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);
1292:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);
1293:   return(0);
1294: }

1296: static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
1297: {
1298:   PetscErrorCode  ierr;
1299:   PetscInt        bs;
1300:   PetscBool       flg;
1301:   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
1302:   PCCompositeType ctype;

1305:   PetscOptionsHead(PetscOptionsObject,"FieldSplit options");
1306:   PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);
1307:   PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);
1308:   if (flg) {
1309:     PCFieldSplitSetBlockSize(pc,bs);
1310:   }
1311:   jac->diag_use_amat = pc->useAmat;
1312:   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);
1313:   jac->offdiag_use_amat = pc->useAmat;
1314:   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);
1315:   PetscOptionsBool("-pc_fieldsplit_detect_saddle_point","Form 2-way split by detecting zero diagonal entries", "PCFieldSplitSetDetectSaddlePoint",jac->detect,&jac->detect,NULL);
1316:   PCFieldSplitSetDetectSaddlePoint(pc,jac->detect); /* Sets split type and Schur PC type */
1317:   PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);
1318:   if (flg) {
1319:     PCFieldSplitSetType(pc,ctype);
1320:   }
1321:   /* Only setup fields once */
1322:   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1323:     /* only allow user to set fields from command line if bs is already known.
1324:        otherwise user can set them in PCFieldSplitSetDefaults() */
1325:     PCFieldSplitSetRuntimeSplits_Private(pc);
1326:     if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
1327:   }
1328:   if (jac->type == PC_COMPOSITE_SCHUR) {
1329:     PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);
1330:     if (flg) {PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");}
1331:     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);
1332:     PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);
1333:     PetscOptionsScalar("-pc_fieldsplit_schur_scale","Scale Schur complement","PCFieldSplitSetSchurScale",jac->schurscale,&jac->schurscale,NULL);
1334:   }
1335:   PetscOptionsTail();
1336:   return(0);
1337: }

1339: /*------------------------------------------------------------------------------------*/

1341: static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1342: {
1343:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1344:   PetscErrorCode    ierr;
1345:   PC_FieldSplitLink ilink,next = jac->head;
1346:   char              prefix[128];
1347:   PetscInt          i;

1350:   if (jac->splitdefined) {
1351:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1352:     return(0);
1353:   }
1354:   for (i=0; i<n; i++) {
1355:     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1356:     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1357:   }
1358:   PetscNew(&ilink);
1359:   if (splitname) {
1360:     PetscStrallocpy(splitname,&ilink->splitname);
1361:   } else {
1362:     PetscMalloc1(3,&ilink->splitname);
1363:     PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);
1364:   }
1365:   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 */
1366:   PetscMalloc1(n,&ilink->fields);
1367:   PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));
1368:   PetscMalloc1(n,&ilink->fields_col);
1369:   PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));

1371:   ilink->nfields = n;
1372:   ilink->next    = NULL;
1373:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1374:   KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1375:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1376:   KSPSetType(ilink->ksp,KSPPREONLY);
1377:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

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

1382:   if (!next) {
1383:     jac->head       = ilink;
1384:     ilink->previous = NULL;
1385:   } else {
1386:     while (next->next) {
1387:       next = next->next;
1388:     }
1389:     next->next      = ilink;
1390:     ilink->previous = next;
1391:   }
1392:   jac->nsplits++;
1393:   return(0);
1394: }

1396: static PetscErrorCode  PCFieldSplitSchurGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1397: {
1398:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1402:   *subksp = NULL;
1403:   if (n) *n = 0;
1404:   if (jac->type == PC_COMPOSITE_SCHUR) {
1405:     PetscInt nn;

1407:     if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitSchurGetSubKSP()");
1408:     if (jac->nsplits != 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unexpected number of splits %D != 2",jac->nsplits);
1409:     nn   = jac->nsplits + (jac->kspupper != jac->head->ksp ? 1 : 0);
1410:     PetscMalloc1(nn,subksp);
1411:     (*subksp)[0] = jac->head->ksp;
1412:     (*subksp)[1] = jac->kspschur;
1413:     if (jac->kspupper != jac->head->ksp) (*subksp)[2] = jac->kspupper;
1414:     if (n) *n = nn;
1415:   }
1416:   return(0);
1417: }

1419: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1420: {
1421:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1425:   if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()");
1426:   PetscMalloc1(jac->nsplits,subksp);
1427:   MatSchurComplementGetKSP(jac->schur,*subksp);

1429:   (*subksp)[1] = jac->kspschur;
1430:   if (n) *n = jac->nsplits;
1431:   return(0);
1432: }

1434: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1435: {
1436:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1437:   PetscErrorCode    ierr;
1438:   PetscInt          cnt   = 0;
1439:   PC_FieldSplitLink ilink = jac->head;

1442:   PetscMalloc1(jac->nsplits,subksp);
1443:   while (ilink) {
1444:     (*subksp)[cnt++] = ilink->ksp;
1445:     ilink            = ilink->next;
1446:   }
1447:   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);
1448:   if (n) *n = jac->nsplits;
1449:   return(0);
1450: }

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

1455:     Input Parameters:
1456: +   pc  - the preconditioner context
1457: +   is - the index set that defines the indices to which the fieldsplit is to be restricted

1459:     Level: advanced

1461: @*/
1462: PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
1463: {

1469:   PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));
1470:   return(0);
1471: }


1474: static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
1475: {
1476:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1477:   PetscErrorCode    ierr;
1478:   PC_FieldSplitLink ilink = jac->head, next;
1479:   PetscInt          localsize,size,sizez,i;
1480:   const PetscInt    *ind, *indz;
1481:   PetscInt          *indc, *indcz;
1482:   PetscBool         flg;

1485:   ISGetLocalSize(isy,&localsize);
1486:   MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));
1487:   size -= localsize;
1488:   while(ilink) {
1489:     IS isrl,isr;
1490:     PC subpc;
1491:     ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);
1492:     ISGetLocalSize(isrl,&localsize);
1493:     PetscMalloc1(localsize,&indc);
1494:     ISGetIndices(isrl,&ind);
1495:     PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));
1496:     ISRestoreIndices(isrl,&ind);
1497:     ISDestroy(&isrl);
1498:     for (i=0; i<localsize; i++) *(indc+i) += size;
1499:     ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);
1500:     PetscObjectReference((PetscObject)isr);
1501:     ISDestroy(&ilink->is);
1502:     ilink->is     = isr;
1503:     PetscObjectReference((PetscObject)isr);
1504:     ISDestroy(&ilink->is_col);
1505:     ilink->is_col = isr;
1506:     ISDestroy(&isr);
1507:     KSPGetPC(ilink->ksp, &subpc);
1508:     PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);
1509:     if(flg) {
1510:       IS iszl,isz;
1511:       MPI_Comm comm;
1512:       ISGetLocalSize(ilink->is,&localsize);
1513:       comm   = PetscObjectComm((PetscObject)ilink->is);
1514:       ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);
1515:       MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);
1516:       sizez -= localsize;
1517:       ISGetLocalSize(iszl,&localsize);
1518:       PetscMalloc1(localsize,&indcz);
1519:       ISGetIndices(iszl,&indz);
1520:       PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));
1521:       ISRestoreIndices(iszl,&indz);
1522:       ISDestroy(&iszl);
1523:       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
1524:       ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);
1525:       PCFieldSplitRestrictIS(subpc,isz);
1526:       ISDestroy(&isz);
1527:     }
1528:     next = ilink->next;
1529:     ilink = next;
1530:   }
1531:   jac->isrestrict = PETSC_TRUE;
1532:   return(0);
1533: }

1535: static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1536: {
1537:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1538:   PetscErrorCode    ierr;
1539:   PC_FieldSplitLink ilink, next = jac->head;
1540:   char              prefix[128];

1543:   if (jac->splitdefined) {
1544:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1545:     return(0);
1546:   }
1547:   PetscNew(&ilink);
1548:   if (splitname) {
1549:     PetscStrallocpy(splitname,&ilink->splitname);
1550:   } else {
1551:     PetscMalloc1(8,&ilink->splitname);
1552:     PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);
1553:   }
1554:   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 */
1555:   PetscObjectReference((PetscObject)is);
1556:   ISDestroy(&ilink->is);
1557:   ilink->is     = is;
1558:   PetscObjectReference((PetscObject)is);
1559:   ISDestroy(&ilink->is_col);
1560:   ilink->is_col = is;
1561:   ilink->next   = NULL;
1562:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1563:   KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1564:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1565:   KSPSetType(ilink->ksp,KSPPREONLY);
1566:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

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

1571:   if (!next) {
1572:     jac->head       = ilink;
1573:     ilink->previous = NULL;
1574:   } else {
1575:     while (next->next) {
1576:       next = next->next;
1577:     }
1578:     next->next      = ilink;
1579:     ilink->previous = next;
1580:   }
1581:   jac->nsplits++;
1582:   return(0);
1583: }

1585: /*@C
1586:     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner

1588:     Logically Collective on PC

1590:     Input Parameters:
1591: +   pc  - the preconditioner context
1592: .   splitname - name of this split, if NULL the number of the split is used
1593: .   n - the number of fields in this split
1594: -   fields - the fields in this split

1596:     Level: intermediate

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

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

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

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

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

1615: @*/
1616: PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1617: {

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

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

1632:     Logically Collective on PC

1634:     Input Parameters:
1635: +   pc  - the preconditioner object
1636: -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from

1638:     Options Database:
1639: .     -pc_fieldsplit_diag_use_amat

1641:     Level: intermediate

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

1645: @*/
1646: PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1647: {
1648:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1649:   PetscBool      isfs;

1654:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1655:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1656:   jac->diag_use_amat = flg;
1657:   return(0);
1658: }

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

1663:     Logically Collective on PC

1665:     Input Parameters:
1666: .   pc  - the preconditioner object

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


1672:     Level: intermediate

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

1676: @*/
1677: PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1678: {
1679:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1680:   PetscBool      isfs;

1686:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1687:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1688:   *flg = jac->diag_use_amat;
1689:   return(0);
1690: }

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

1695:     Logically Collective on PC

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

1701:     Options Database:
1702: .     -pc_fieldsplit_off_diag_use_amat

1704:     Level: intermediate

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

1708: @*/
1709: PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1710: {
1711:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1712:   PetscBool      isfs;

1717:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1718:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1719:   jac->offdiag_use_amat = flg;
1720:   return(0);
1721: }

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

1726:     Logically Collective on PC

1728:     Input Parameters:
1729: .   pc  - the preconditioner object

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


1735:     Level: intermediate

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

1739: @*/
1740: PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1741: {
1742:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1743:   PetscBool      isfs;

1749:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1750:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1751:   *flg = jac->offdiag_use_amat;
1752:   return(0);
1753: }



1757: /*@C
1758:     PCFieldSplitSetIS - Sets the exact elements for field

1760:     Logically Collective on PC

1762:     Input Parameters:
1763: +   pc  - the preconditioner context
1764: .   splitname - name of this split, if NULL the number of the split is used
1765: -   is - the index set that defines the vector elements in this field


1768:     Notes:
1769:     Use PCFieldSplitSetFields(), for fields defined by strided types.

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

1774:     Level: intermediate

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

1778: @*/
1779: PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1780: {

1787:   PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));
1788:   return(0);
1789: }

1791: /*@C
1792:     PCFieldSplitGetIS - Retrieves the elements for a field as an IS

1794:     Logically Collective on PC

1796:     Input Parameters:
1797: +   pc  - the preconditioner context
1798: -   splitname - name of this split

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

1803:     Level: intermediate

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

1807: @*/
1808: PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
1809: {

1816:   {
1817:     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
1818:     PC_FieldSplitLink ilink = jac->head;
1819:     PetscBool         found;

1821:     *is = NULL;
1822:     while (ilink) {
1823:       PetscStrcmp(ilink->splitname, splitname, &found);
1824:       if (found) {
1825:         *is = ilink->is;
1826:         break;
1827:       }
1828:       ilink = ilink->next;
1829:     }
1830:   }
1831:   return(0);
1832: }

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

1838:     Logically Collective on PC

1840:     Input Parameters:
1841: +   pc  - the preconditioner context
1842: -   bs - the block size

1844:     Level: intermediate

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

1848: @*/
1849: PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1850: {

1856:   PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));
1857:   return(0);
1858: }

1860: /*@C
1861:    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits

1863:    Collective on KSP

1865:    Input Parameter:
1866: .  pc - the preconditioner context

1868:    Output Parameters:
1869: +  n - the number of splits
1870: -  subksp - the array of KSP contexts

1872:    Note:
1873:    After PCFieldSplitGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
1874:    (not the KSP just the array that contains them).

1876:    You must call PCSetUp() before calling PCFieldSplitGetSubKSP().

1878:    If the fieldsplit is of type PC_COMPOSITE_SCHUR, it returns the KSP object used inside the
1879:    Schur complement and the KSP object used to iterate over the Schur complement.
1880:    To access all the KSP objects used in PC_COMPOSITE_SCHUR, use PCFieldSplitSchurGetSubKSP()

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


1887:    Level: advanced

1889: .seealso: PCFIELDSPLIT
1890: @*/
1891: PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1892: {

1898:   PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
1899:   return(0);
1900: }

1902: /*@C
1903:    PCFieldSplitSchurGetSubKSP - Gets the KSP contexts used inside the Schur complement based PCFIELDSPLIT

1905:    Collective on KSP

1907:    Input Parameter:
1908: .  pc - the preconditioner context

1910:    Output Parameters:
1911: +  n - the number of splits
1912: -  subksp - the array of KSP contexts

1914:    Note:
1915:    After PCFieldSplitSchurGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
1916:    (not the KSP just the array that contains them).

1918:    You must call PCSetUp() before calling PCFieldSplitSchurGetSubKSP().

1920:    If the fieldsplit type is of type PC_COMPOSITE_SCHUR, it returns (in order)
1921:    - the KSP used for the (1,1) block
1922:    - the KSP used for the Schur complement (not the one used for the interior Schur solver)
1923:    - the KSP used for the (1,1) block in the upper triangular factor (if different from that of the (1,1) block).

1925:    It returns a null array if the fieldsplit is not of type PC_COMPOSITE_SCHUR; in this case, you should use PCFieldSplitGetSubKSP().

1927:    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1928:       You can call PCFieldSplitSchurGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
1929:       KSP array must be.

1931:    Level: advanced

1933: .seealso: PCFIELDSPLIT
1934: @*/
1935: PetscErrorCode  PCFieldSplitSchurGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1936: {

1942:   PetscUseMethod(pc,"PCFieldSplitSchurGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
1943:   return(0);
1944: }

1946: /*@
1947:     PCFieldSplitSetSchurPre -  Indicates what operator is used to construct the preconditioner for the Schur complement.
1948:       A11 matrix. Otherwise no preconditioner is used.

1950:     Collective on PC

1952:     Input Parameters:
1953: +   pc      - the preconditioner context
1954: .   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_SCHUR_PRE_USER
1955:               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
1956: -   userpre - matrix to use for preconditioning, or NULL

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

1961:     Notes:
1962: $    If ptype is
1963: $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
1964: $             matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
1965: $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
1966: $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
1967: $             preconditioner
1968: $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
1969: $             to this function).
1970: $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
1971: $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
1972: $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
1973: $        full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PCFIELDSPLIT (this is expensive)
1974: $             useful mostly as a test that the Schur complement approach can work for your problem

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

1980:     Level: intermediate

1982: .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
1983:           MatSchurComplementSetAinvType(), PCLSC

1985: @*/
1986: PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1987: {

1992:   PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));
1993:   return(0);
1994: }

1996: PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */

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

2002:     Logically Collective on PC

2004:     Input Parameters:
2005: .   pc      - the preconditioner context

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

2011:     Level: intermediate

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

2015: @*/
2016: PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
2017: {

2022:   PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));
2023:   return(0);
2024: }

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

2029:     Not collective

2031:     Input Parameter:
2032: .   pc      - the preconditioner context

2034:     Output Parameter:
2035: .   S       - the Schur complement matrix

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

2040:     Level: advanced

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

2044: @*/
2045: PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
2046: {
2048:   const char*    t;
2049:   PetscBool      isfs;
2050:   PC_FieldSplit  *jac;

2054:   PetscObjectGetType((PetscObject)pc,&t);
2055:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
2056:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2057:   jac = (PC_FieldSplit*)pc->data;
2058:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
2059:   if (S) *S = jac->schur;
2060:   return(0);
2061: }

2063: /*@
2064:     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC

2066:     Not collective

2068:     Input Parameters:
2069: +   pc      - the preconditioner context
2070: .   S       - the Schur complement matrix

2072:     Level: advanced

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

2076: @*/
2077: PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
2078: {
2080:   const char*    t;
2081:   PetscBool      isfs;
2082:   PC_FieldSplit  *jac;

2086:   PetscObjectGetType((PetscObject)pc,&t);
2087:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
2088:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2089:   jac = (PC_FieldSplit*)pc->data;
2090:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
2091:   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
2092:   return(0);
2093: }


2096: static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2097: {
2098:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2102:   jac->schurpre = ptype;
2103:   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
2104:     MatDestroy(&jac->schur_user);
2105:     jac->schur_user = pre;
2106:     PetscObjectReference((PetscObject)jac->schur_user);
2107:   }
2108:   return(0);
2109: }

2111: static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
2112: {
2113:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2116:   *ptype = jac->schurpre;
2117:   *pre   = jac->schur_user;
2118:   return(0);
2119: }

2121: /*@
2122:     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain in the preconditioner

2124:     Collective on PC

2126:     Input Parameters:
2127: +   pc  - the preconditioner context
2128: -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default

2130:     Options Database:
2131: .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full


2134:     Level: intermediate

2136:     Notes:
2137:     The FULL factorization is

2139: $   (A   B)  = (1       0) (A   0) (1  Ainv*B)  = L D U
2140: $   (C   E)    (C*Ainv  1) (0   S) (0     1  )

2142:     where S = E - 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,
2143:     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). Sign flipping of S can be turned off with PCFieldSplitSetSchurScale().

2145: $    If A and S are solved exactly
2146: $      *) FULL factorization is a direct solver.
2147: $      *) The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial of degree 2, so KSPGMRES converges in 2 iterations.
2148: $      *) With DIAG, the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.

2150:     If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
2151:     application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice.

2153:     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES.

2155:     Note that a flexible method like KSPFGMRES or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used to solve with A or S).

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

2161: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCFieldSplitSetSchurScale()
2162: @*/
2163: PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2164: {

2169:   PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));
2170:   return(0);
2171: }

2173: static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2174: {
2175:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2178:   jac->schurfactorization = ftype;
2179:   return(0);
2180: }

2182: /*@
2183:     PCFieldSplitSetSchurScale -  Controls the sign flip of S for PC_FIELDSPLIT_SCHUR_FACT_DIAG.

2185:     Collective on PC

2187:     Input Parameters:
2188: +   pc    - the preconditioner context
2189: -   scale - scaling factor for the Schur complement

2191:     Options Database:
2192: .     -pc_fieldsplit_schur_scale - default is -1.0

2194:     Level: intermediate

2196: .seealso: PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurFactType, PCFieldSplitSetSchurScale()
2197: @*/
2198: PetscErrorCode PCFieldSplitSetSchurScale(PC pc,PetscScalar scale)
2199: {

2205:   PetscTryMethod(pc,"PCFieldSplitSetSchurScale_C",(PC,PetscScalar),(pc,scale));
2206:   return(0);
2207: }

2209: static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc,PetscScalar scale)
2210: {
2211:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2214:   jac->schurscale = scale;
2215:   return(0);
2216: }

2218: /*@C
2219:    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement

2221:    Collective on KSP

2223:    Input Parameter:
2224: .  pc - the preconditioner context

2226:    Output Parameters:
2227: +  A00 - the (0,0) block
2228: .  A01 - the (0,1) block
2229: .  A10 - the (1,0) block
2230: -  A11 - the (1,1) block

2232:    Level: advanced

2234: .seealso: PCFIELDSPLIT
2235: @*/
2236: PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
2237: {
2238:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2242:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2243:   if (A00) *A00 = jac->pmat[0];
2244:   if (A01) *A01 = jac->B;
2245:   if (A10) *A10 = jac->C;
2246:   if (A11) *A11 = jac->pmat[1];
2247:   return(0);
2248: }

2250: static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
2251: {
2252:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2256:   jac->type = type;
2257:   if (type == PC_COMPOSITE_SCHUR) {
2258:     pc->ops->apply = PCApply_FieldSplit_Schur;
2259:     pc->ops->view  = PCView_FieldSplit_Schur;

2261:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);
2262:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);
2263:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);
2264:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);
2265:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",PCFieldSplitSetSchurScale_FieldSplit);

2267:   } else {
2268:     pc->ops->apply = PCApply_FieldSplit;
2269:     pc->ops->view  = PCView_FieldSplit;

2271:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2272:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);
2273:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);
2274:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);
2275:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",0);
2276:   }
2277:   return(0);
2278: }

2280: static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
2281: {
2282:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2285:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2286:   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);
2287:   jac->bs = bs;
2288:   return(0);
2289: }

2291: /*@
2292:    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.

2294:    Collective on PC

2296:    Input Parameter:
2297: .  pc - the preconditioner context
2298: .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

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

2303:    Level: Intermediate

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

2307: .seealso: PCCompositeSetType()

2309: @*/
2310: PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
2311: {

2316:   PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));
2317:   return(0);
2318: }

2320: /*@
2321:   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.

2323:   Not collective

2325:   Input Parameter:
2326: . pc - the preconditioner context

2328:   Output Parameter:
2329: . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

2331:   Level: Intermediate

2333: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2334: .seealso: PCCompositeSetType()
2335: @*/
2336: PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2337: {
2338:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2343:   *type = jac->type;
2344:   return(0);
2345: }

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

2350:    Logically Collective

2352:    Input Parameters:
2353: +  pc   - the preconditioner context
2354: -  flg  - boolean indicating whether to use field splits defined by the DM

2356:    Options Database Key:
2357: .  -pc_fieldsplit_dm_splits

2359:    Level: Intermediate

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

2363: .seealso: PCFieldSplitGetDMSplits()

2365: @*/
2366: PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
2367: {
2368:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2369:   PetscBool      isfs;

2375:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2376:   if (isfs) {
2377:     jac->dm_splits = flg;
2378:   }
2379:   return(0);
2380: }


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

2386:    Logically Collective

2388:    Input Parameter:
2389: .  pc   - the preconditioner context

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

2394:    Level: Intermediate

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

2398: .seealso: PCFieldSplitSetDMSplits()

2400: @*/
2401: PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
2402: {
2403:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2404:   PetscBool      isfs;

2410:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2411:   if (isfs) {
2412:     if(flg) *flg = jac->dm_splits;
2413:   }
2414:   return(0);
2415: }

2417: /*@
2418:    PCFieldSplitGetDetectSaddlePoint - Returns flag indicating whether PCFieldSplit will attempt to automatically determine fields based on zero diagonal entries.

2420:    Logically Collective

2422:    Input Parameter:
2423: .  pc   - the preconditioner context

2425:    Output Parameter:
2426: .  flg  - boolean indicating whether to detect fields or not

2428:    Level: Intermediate

2430: .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint()

2432: @*/
2433: PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC pc,PetscBool *flg)
2434: {
2435:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2438:   *flg = jac->detect;
2439:   return(0);
2440: }

2442: /*@
2443:    PCFieldSplitSetDetectSaddlePoint - Sets flag indicating whether PCFieldSplit will attempt to automatically determine fields based on zero diagonal entries.

2445:    Logically Collective

2447:    Notes:
2448:    Also sets the split type to PC_COMPOSITE_SCHUR (see PCFieldSplitSetType()) and the Schur preconditioner type to PC_FIELDSPLIT_SCHUR_PRE_SELF (see PCFieldSplitSetSchurPre()).

2450:    Input Parameter:
2451: .  pc   - the preconditioner context

2453:    Output Parameter:
2454: .  flg  - boolean indicating whether to detect fields or not

2456:    Options Database Key:
2457: .  -pc_fieldsplit_detect_saddle_point

2459:    Level: Intermediate

2461: .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint(), PCFieldSplitSetType(), PCFieldSplitSetSchurPre()

2463: @*/
2464: PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc,PetscBool flg)
2465: {
2466:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2470:   jac->detect = flg;
2471:   if (jac->detect) {
2472:     PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
2473:     PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELF,NULL);
2474:   }
2475:   return(0);
2476: }

2478: /* -------------------------------------------------------------------------------------*/
2479: /*MC
2480:    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2481:                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.

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

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

2489:    Level: intermediate

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

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

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

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

2512: $     For the Schur complement preconditioner if J = ( A00 A01 )
2513: $                                                    ( A10 A11 )
2514: $     the preconditioner using full factorization is
2515: $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2516: $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
2517:      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
2518: $              S = A11 - A10 ksp(A00) A01
2519:      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
2520:      in providing the SECOND split or 1 if not give). For PCFieldSplitGetSubKSP() when field number is 0,
2521:      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2522:      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.

2524:      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
2525:      diag gives
2526: $              ( inv(A00)     0   )
2527: $              (   0      -ksp(S) )
2528:      note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. For SPD matrices J, the sign flip
2529:      can be turned off with PCFieldSplitSetSchurScale() or by command line -pc_fieldsplit_schur_scale 1.0. The lower factorization is the inverse of
2530: $              (  A00   0 )
2531: $              (  A10   S )
2532:      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2533: $              ( A00 A01 )
2534: $              (  0   S  )
2535:      where again the inverses of A00 and S are applied using KSPs.

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

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

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

2547:    Concepts: physics based preconditioners, block preconditioners

2549:    There is a nice discussion of block preconditioners in

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

2555:    The Constrained Pressure Preconditioner (CPR) can be implemented using PCCOMPOSITE with PCGALERKIN. CPR first solves an R A P subsystem, updates the
2556:    residual on all variables (PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE)), and then applies a simple ILU like preconditioner on all the variables.

2558: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
2559:            PCFieldSplitGetSubKSP(), PCFieldSplitSchurGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
2560:           MatSchurComplementSetAinvType(), PCFieldSplitSetSchurScale(),
2561:           PCFieldSplitSetDetectSaddlePoint()
2562: M*/

2564: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
2565: {
2567:   PC_FieldSplit  *jac;

2570:   PetscNewLog(pc,&jac);

2572:   jac->bs                 = -1;
2573:   jac->nsplits            = 0;
2574:   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2575:   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2576:   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2577:   jac->schurscale         = -1.0;
2578:   jac->dm_splits          = PETSC_TRUE;
2579:   jac->detect             = PETSC_FALSE;

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

2583:   pc->ops->apply           = PCApply_FieldSplit;
2584:   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
2585:   pc->ops->setup           = PCSetUp_FieldSplit;
2586:   pc->ops->reset           = PCReset_FieldSplit;
2587:   pc->ops->destroy         = PCDestroy_FieldSplit;
2588:   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
2589:   pc->ops->view            = PCView_FieldSplit;
2590:   pc->ops->applyrichardson = 0;

2592:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurGetSubKSP_C",PCFieldSplitSchurGetSubKSP_FieldSplit);
2593:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2594:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);
2595:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);
2596:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);
2597:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);
2598:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);
2599:   return(0);
2600: }