Actual source code: fieldsplit.c

petsc-master 2018-04-26
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:     PetscViewerASCIIPushTab(viewer);
110:     for (i=0; i<jac->nsplits; i++) {
111:       if (ilink->fields) {
112:         PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
113:         PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
114:         for (j=0; j<ilink->nfields; j++) {
115:           if (j > 0) {
116:             PetscViewerASCIIPrintf(viewer,",");
117:           }
118:           PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
119:         }
120:         PetscViewerASCIIPrintf(viewer,"\n");
121:         PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
122:       } else {
123:         PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);
124:       }
125:       KSPView(ilink->ksp,viewer);
126:       ilink = ilink->next;
127:     }
128:     PetscViewerASCIIPopTab(viewer);
129:   }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

472:   PCFieldSplitSetDefaults(pc);
473:   nsplit = jac->nsplits;
474:   ilink  = jac->head;

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

480:     jac->issetup = PETSC_TRUE;

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

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

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

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

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

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

588:   /* Check for null space attached to IS */
589:   ilink = jac->head;
590:   for (i=0; i<nsplit; i++) {
591:     MatNullSpace sp;

593:     PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);
594:     if (sp) {
595:       MatSetNullSpace(jac->mat[i], sp);
596:     }
597:     ilink = ilink->next;
598:   }

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

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

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

653:     /* If pc->mat is SPD, don't scale by -1 the Schur complement */
654:     if (jac->schurscale == (PetscScalar)-1.0) {
655:       MatGetOption(pc->pmat,MAT_SPD,&isspd);
656:       jac->schurscale = (isspd == PETSC_TRUE) ? 1.0 : -1.0;
657:     }

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

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

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

702:       /* extract the A01 and A10 matrices */
703:       ilink = jac->head;
704:       ISComplement(ilink->is_col,rstart,rend,&ccis);
705:       if (jac->offdiag_use_amat) {
706:         MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
707:       } else {
708:         MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
709:       }
710:       ISDestroy(&ccis);
711:       ilink = ilink->next;
712:       ISComplement(ilink->is_col,rstart,rend,&ccis);
713:       if (jac->offdiag_use_amat) {
714:         MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
715:       } else {
716:         MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
717:       }
718:       ISDestroy(&ccis);

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

733:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);
734:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
735:       if (flg) {
736:         DM  dmInner;
737:         KSP kspInner;

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

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

763:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);
764:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
765:       if (flg) {
766:         DM dmInner;

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

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

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

837:   jac->suboptionsset = PETSC_TRUE;
838:   return(0);
839: }

841: #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
842:   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
843:    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
844:    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
845:    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
846:    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
847:    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
848:    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))

850: static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
851: {
852:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
853:   PetscErrorCode    ierr;
854:   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
855:   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;

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

926:     PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
927:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
928:     PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
929:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
930:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);

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

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

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

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

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

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

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

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

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

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

1213: static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1214: {
1215:   PetscErrorCode    ierr;

1218:   PCReset_FieldSplit(pc);
1219:   PetscFree(pc->data);
1220:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);
1221:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);
1222:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);
1223:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);
1224:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);
1225:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);
1226:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);
1227:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);
1228:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);
1229:   return(0);
1230: }

1232: static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
1233: {
1234:   PetscErrorCode  ierr;
1235:   PetscInt        bs;
1236:   PetscBool       flg;
1237:   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
1238:   PCCompositeType ctype;

1241:   PetscOptionsHead(PetscOptionsObject,"FieldSplit options");
1242:   PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);
1243:   PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);
1244:   if (flg) {
1245:     PCFieldSplitSetBlockSize(pc,bs);
1246:   }
1247:   jac->diag_use_amat = pc->useAmat;
1248:   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);
1249:   jac->offdiag_use_amat = pc->useAmat;
1250:   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);
1251:   PetscOptionsBool("-pc_fieldsplit_detect_saddle_point","Form 2-way split by detecting zero diagonal entries", "PCFieldSplitSetDetectSaddlePoint",jac->detect,&jac->detect,NULL);
1252:   PCFieldSplitSetDetectSaddlePoint(pc,jac->detect); /* Sets split type and Schur PC type */
1253:   PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);
1254:   if (flg) {
1255:     PCFieldSplitSetType(pc,ctype);
1256:   }
1257:   /* Only setup fields once */
1258:   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1259:     /* only allow user to set fields from command line if bs is already known.
1260:        otherwise user can set them in PCFieldSplitSetDefaults() */
1261:     PCFieldSplitSetRuntimeSplits_Private(pc);
1262:     if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
1263:   }
1264:   if (jac->type == PC_COMPOSITE_SCHUR) {
1265:     PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);
1266:     if (flg) {PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");}
1267:     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);
1268:     PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);
1269:     PetscOptionsScalar("-pc_fieldsplit_schur_scale","Scale Schur complement","PCFieldSplitSetSchurScale",jac->schurscale,&jac->schurscale,NULL);
1270:   }
1271:   PetscOptionsTail();
1272:   return(0);
1273: }

1275: /*------------------------------------------------------------------------------------*/

1277: static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1278: {
1279:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1280:   PetscErrorCode    ierr;
1281:   PC_FieldSplitLink ilink,next = jac->head;
1282:   char              prefix[128];
1283:   PetscInt          i;

1286:   if (jac->splitdefined) {
1287:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1288:     return(0);
1289:   }
1290:   for (i=0; i<n; i++) {
1291:     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1292:     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1293:   }
1294:   PetscNew(&ilink);
1295:   if (splitname) {
1296:     PetscStrallocpy(splitname,&ilink->splitname);
1297:   } else {
1298:     PetscMalloc1(3,&ilink->splitname);
1299:     PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);
1300:   }
1301:   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 */
1302:   PetscMalloc1(n,&ilink->fields);
1303:   PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));
1304:   PetscMalloc1(n,&ilink->fields_col);
1305:   PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));

1307:   ilink->nfields = n;
1308:   ilink->next    = NULL;
1309:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1310:   KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1311:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1312:   KSPSetType(ilink->ksp,KSPPREONLY);
1313:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

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

1318:   if (!next) {
1319:     jac->head       = ilink;
1320:     ilink->previous = NULL;
1321:   } else {
1322:     while (next->next) {
1323:       next = next->next;
1324:     }
1325:     next->next      = ilink;
1326:     ilink->previous = next;
1327:   }
1328:   jac->nsplits++;
1329:   return(0);
1330: }

1332: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1333: {
1334:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

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

1342:   (*subksp)[1] = jac->kspschur;
1343:   if (n) *n = jac->nsplits;
1344:   return(0);
1345: }

1347: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1348: {
1349:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1350:   PetscErrorCode    ierr;
1351:   PetscInt          cnt   = 0;
1352:   PC_FieldSplitLink ilink = jac->head;

1355:   PetscMalloc1(jac->nsplits,subksp);
1356:   while (ilink) {
1357:     (*subksp)[cnt++] = ilink->ksp;
1358:     ilink            = ilink->next;
1359:   }
1360:   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);
1361:   if (n) *n = jac->nsplits;
1362:   return(0);
1363: }

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

1368:     Input Parameters:
1369: +   pc  - the preconditioner context
1370: +   is - the index set that defines the indices to which the fieldsplit is to be restricted

1372:     Level: advanced

1374: @*/
1375: PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
1376: {

1382:   PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));
1383:   return(0);
1384: }


1387: static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
1388: {
1389:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1390:   PetscErrorCode    ierr;
1391:   PC_FieldSplitLink ilink = jac->head, next;
1392:   PetscInt          localsize,size,sizez,i;
1393:   const PetscInt    *ind, *indz;
1394:   PetscInt          *indc, *indcz;
1395:   PetscBool         flg;

1398:   ISGetLocalSize(isy,&localsize);
1399:   MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));
1400:   size -= localsize;
1401:   while(ilink) {
1402:     IS isrl,isr;
1403:     PC subpc;
1404:     ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);
1405:     ISGetLocalSize(isrl,&localsize);
1406:     PetscMalloc1(localsize,&indc);
1407:     ISGetIndices(isrl,&ind);
1408:     PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));
1409:     ISRestoreIndices(isrl,&ind);
1410:     ISDestroy(&isrl);
1411:     for (i=0; i<localsize; i++) *(indc+i) += size;
1412:     ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);
1413:     PetscObjectReference((PetscObject)isr);
1414:     ISDestroy(&ilink->is);
1415:     ilink->is     = isr;
1416:     PetscObjectReference((PetscObject)isr);
1417:     ISDestroy(&ilink->is_col);
1418:     ilink->is_col = isr;
1419:     ISDestroy(&isr);
1420:     KSPGetPC(ilink->ksp, &subpc);
1421:     PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);
1422:     if(flg) {
1423:       IS iszl,isz;
1424:       MPI_Comm comm;
1425:       ISGetLocalSize(ilink->is,&localsize);
1426:       comm   = PetscObjectComm((PetscObject)ilink->is);
1427:       ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);
1428:       MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);
1429:       sizez -= localsize;
1430:       ISGetLocalSize(iszl,&localsize);
1431:       PetscMalloc1(localsize,&indcz);
1432:       ISGetIndices(iszl,&indz);
1433:       PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));
1434:       ISRestoreIndices(iszl,&indz);
1435:       ISDestroy(&iszl);
1436:       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
1437:       ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);
1438:       PCFieldSplitRestrictIS(subpc,isz);
1439:       ISDestroy(&isz);
1440:     }
1441:     next = ilink->next;
1442:     ilink = next;
1443:   }
1444:   jac->isrestrict = PETSC_TRUE;
1445:   return(0);
1446: }

1448: static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1449: {
1450:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1451:   PetscErrorCode    ierr;
1452:   PC_FieldSplitLink ilink, next = jac->head;
1453:   char              prefix[128];

1456:   if (jac->splitdefined) {
1457:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1458:     return(0);
1459:   }
1460:   PetscNew(&ilink);
1461:   if (splitname) {
1462:     PetscStrallocpy(splitname,&ilink->splitname);
1463:   } else {
1464:     PetscMalloc1(8,&ilink->splitname);
1465:     PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);
1466:   }
1467:   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 */
1468:   PetscObjectReference((PetscObject)is);
1469:   ISDestroy(&ilink->is);
1470:   ilink->is     = is;
1471:   PetscObjectReference((PetscObject)is);
1472:   ISDestroy(&ilink->is_col);
1473:   ilink->is_col = is;
1474:   ilink->next   = NULL;
1475:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1476:   KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1477:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1478:   KSPSetType(ilink->ksp,KSPPREONLY);
1479:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

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

1484:   if (!next) {
1485:     jac->head       = ilink;
1486:     ilink->previous = NULL;
1487:   } else {
1488:     while (next->next) {
1489:       next = next->next;
1490:     }
1491:     next->next      = ilink;
1492:     ilink->previous = next;
1493:   }
1494:   jac->nsplits++;
1495:   return(0);
1496: }

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

1501:     Logically Collective on PC

1503:     Input Parameters:
1504: +   pc  - the preconditioner context
1505: .   splitname - name of this split, if NULL the number of the split is used
1506: .   n - the number of fields in this split
1507: -   fields - the fields in this split

1509:     Level: intermediate

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

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

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

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

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

1528: @*/
1529: PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1530: {

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

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

1545:     Logically Collective on PC

1547:     Input Parameters:
1548: +   pc  - the preconditioner object
1549: -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from

1551:     Options Database:
1552: .     -pc_fieldsplit_diag_use_amat

1554:     Level: intermediate

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

1558: @*/
1559: PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1560: {
1561:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1562:   PetscBool      isfs;

1567:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1568:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1569:   jac->diag_use_amat = flg;
1570:   return(0);
1571: }

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

1576:     Logically Collective on PC

1578:     Input Parameters:
1579: .   pc  - the preconditioner object

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


1585:     Level: intermediate

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

1589: @*/
1590: PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1591: {
1592:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1593:   PetscBool      isfs;

1599:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1600:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1601:   *flg = jac->diag_use_amat;
1602:   return(0);
1603: }

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

1608:     Logically Collective on PC

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

1614:     Options Database:
1615: .     -pc_fieldsplit_off_diag_use_amat

1617:     Level: intermediate

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

1621: @*/
1622: PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1623: {
1624:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1625:   PetscBool      isfs;

1630:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1631:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1632:   jac->offdiag_use_amat = flg;
1633:   return(0);
1634: }

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

1639:     Logically Collective on PC

1641:     Input Parameters:
1642: .   pc  - the preconditioner object

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


1648:     Level: intermediate

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

1652: @*/
1653: PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1654: {
1655:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1656:   PetscBool      isfs;

1662:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1663:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1664:   *flg = jac->offdiag_use_amat;
1665:   return(0);
1666: }



1670: /*@C
1671:     PCFieldSplitSetIS - Sets the exact elements for field

1673:     Logically Collective on PC

1675:     Input Parameters:
1676: +   pc  - the preconditioner context
1677: .   splitname - name of this split, if NULL the number of the split is used
1678: -   is - the index set that defines the vector elements in this field


1681:     Notes:
1682:     Use PCFieldSplitSetFields(), for fields defined by strided types.

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

1687:     Level: intermediate

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

1691: @*/
1692: PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1693: {

1700:   PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));
1701:   return(0);
1702: }

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

1707:     Logically Collective on PC

1709:     Input Parameters:
1710: +   pc  - the preconditioner context
1711: -   splitname - name of this split

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

1716:     Level: intermediate

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

1720: @*/
1721: PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
1722: {

1729:   {
1730:     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
1731:     PC_FieldSplitLink ilink = jac->head;
1732:     PetscBool         found;

1734:     *is = NULL;
1735:     while (ilink) {
1736:       PetscStrcmp(ilink->splitname, splitname, &found);
1737:       if (found) {
1738:         *is = ilink->is;
1739:         break;
1740:       }
1741:       ilink = ilink->next;
1742:     }
1743:   }
1744:   return(0);
1745: }

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

1751:     Logically Collective on PC

1753:     Input Parameters:
1754: +   pc  - the preconditioner context
1755: -   bs - the block size

1757:     Level: intermediate

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

1761: @*/
1762: PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1763: {

1769:   PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));
1770:   return(0);
1771: }

1773: /*@C
1774:    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits

1776:    Collective on KSP

1778:    Input Parameter:
1779: .  pc - the preconditioner context

1781:    Output Parameters:
1782: +  n - the number of splits
1783: -  subksp - the array of KSP contexts

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

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

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


1796:    Level: advanced

1798: .seealso: PCFIELDSPLIT
1799: @*/
1800: PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1801: {

1807:   PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
1808:   return(0);
1809: }

1811: /*@
1812:     PCFieldSplitSetSchurPre -  Indicates what operator is used to construct the preconditioner for the Schur complement.
1813:       A11 matrix. Otherwise no preconditioner is used.

1815:     Collective on PC

1817:     Input Parameters:
1818: +   pc      - the preconditioner context
1819: .   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 
1820:               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
1821: -   userpre - matrix to use for preconditioning, or NULL

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

1826:     Notes:
1827: $    If ptype is
1828: $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
1829: $             matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
1830: $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
1831: $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
1832: $             preconditioner
1833: $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
1834: $             to this function).
1835: $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
1836: $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
1837: $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
1838: $        full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PFIELDSPLIT (this is expensive)
1839: $             useful mostly as a test that the Schur complement approach can work for your problem

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

1845:     Level: intermediate

1847: .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
1848:           MatSchurComplementSetAinvType(), PCLSC

1850: @*/
1851: PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1852: {

1857:   PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));
1858:   return(0);
1859: }
1860: PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */

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

1866:     Logically Collective on PC

1868:     Input Parameters:
1869: .   pc      - the preconditioner context

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

1875:     Level: intermediate

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

1879: @*/
1880: PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
1881: {

1886:   PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));
1887:   return(0);
1888: }

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

1893:     Not collective

1895:     Input Parameter:
1896: .   pc      - the preconditioner context

1898:     Output Parameter:
1899: .   S       - the Schur complement matrix

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

1904:     Level: advanced

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

1908: @*/
1909: PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
1910: {
1912:   const char*    t;
1913:   PetscBool      isfs;
1914:   PC_FieldSplit  *jac;

1918:   PetscObjectGetType((PetscObject)pc,&t);
1919:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
1920:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1921:   jac = (PC_FieldSplit*)pc->data;
1922:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1923:   if (S) *S = jac->schur;
1924:   return(0);
1925: }

1927: /*@
1928:     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC

1930:     Not collective

1932:     Input Parameters:
1933: +   pc      - the preconditioner context
1934: .   S       - the Schur complement matrix

1936:     Level: advanced

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

1940: @*/
1941: PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
1942: {
1944:   const char*    t;
1945:   PetscBool      isfs;
1946:   PC_FieldSplit  *jac;

1950:   PetscObjectGetType((PetscObject)pc,&t);
1951:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
1952:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1953:   jac = (PC_FieldSplit*)pc->data;
1954:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1955:   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
1956:   return(0);
1957: }


1960: static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1961: {
1962:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1966:   jac->schurpre = ptype;
1967:   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
1968:     MatDestroy(&jac->schur_user);
1969:     jac->schur_user = pre;
1970:     PetscObjectReference((PetscObject)jac->schur_user);
1971:   }
1972:   return(0);
1973: }

1975: static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
1976: {
1977:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1980:   *ptype = jac->schurpre;
1981:   *pre   = jac->schur_user;
1982:   return(0);
1983: }

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

1988:     Collective on PC

1990:     Input Parameters:
1991: +   pc  - the preconditioner context
1992: -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default

1994:     Options Database:
1995: .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full


1998:     Level: intermediate

2000:     Notes:
2001:     The FULL factorization is

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

2006:     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,
2007:     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().

2009: $    If A and S are solved exactly
2010: $      *) FULL factorization is a direct solver.
2011: $      *) 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.
2012: $      *) 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.

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

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

2019:     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).

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

2025: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCFieldSplitSetSchurScale()
2026: @*/
2027: PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2028: {

2033:   PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));
2034:   return(0);
2035: }

2037: static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2038: {
2039:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2042:   jac->schurfactorization = ftype;
2043:   return(0);
2044: }

2046: /*@
2047:     PCFieldSplitSetSchurScale -  Controls the sign flip of S for PC_FIELDSPLIT_SCHUR_FACT_DIAG.

2049:     Collective on PC

2051:     Input Parameters:
2052: +   pc    - the preconditioner context
2053: -   scale - scaling factor for the Schur complement

2055:     Options Database:
2056: .     -pc_fieldsplit_schur_scale - default is -1.0

2058:     Level: intermediate

2060: .seealso: PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurFactType, PCFieldSplitSetSchurScale()
2061: @*/
2062: PetscErrorCode PCFieldSplitSetSchurScale(PC pc,PetscScalar scale)
2063: {

2069:   PetscTryMethod(pc,"PCFieldSplitSetSchurScale_C",(PC,PetscScalar),(pc,scale));
2070:   return(0);
2071: }

2073: static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc,PetscScalar scale)
2074: {
2075:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2078:   jac->schurscale = scale;
2079:   return(0);
2080: }

2082: /*@C
2083:    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement

2085:    Collective on KSP

2087:    Input Parameter:
2088: .  pc - the preconditioner context

2090:    Output Parameters:
2091: +  A00 - the (0,0) block
2092: .  A01 - the (0,1) block
2093: .  A10 - the (1,0) block
2094: -  A11 - the (1,1) block

2096:    Level: advanced

2098: .seealso: PCFIELDSPLIT
2099: @*/
2100: PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
2101: {
2102:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2106:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2107:   if (A00) *A00 = jac->pmat[0];
2108:   if (A01) *A01 = jac->B;
2109:   if (A10) *A10 = jac->C;
2110:   if (A11) *A11 = jac->pmat[1];
2111:   return(0);
2112: }

2114: static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
2115: {
2116:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2120:   jac->type = type;
2121:   if (type == PC_COMPOSITE_SCHUR) {
2122:     pc->ops->apply = PCApply_FieldSplit_Schur;
2123:     pc->ops->view  = PCView_FieldSplit_Schur;

2125:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);
2126:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);
2127:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);
2128:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);
2129:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",PCFieldSplitSetSchurScale_FieldSplit);

2131:   } else {
2132:     pc->ops->apply = PCApply_FieldSplit;
2133:     pc->ops->view  = PCView_FieldSplit;

2135:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2136:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);
2137:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);
2138:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);
2139:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",0);
2140:   }
2141:   return(0);
2142: }

2144: static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
2145: {
2146:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2149:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2150:   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);
2151:   jac->bs = bs;
2152:   return(0);
2153: }

2155: /*@
2156:    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.

2158:    Collective on PC

2160:    Input Parameter:
2161: .  pc - the preconditioner context
2162: .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

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

2167:    Level: Intermediate

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

2171: .seealso: PCCompositeSetType()

2173: @*/
2174: PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
2175: {

2180:   PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));
2181:   return(0);
2182: }

2184: /*@
2185:   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.

2187:   Not collective

2189:   Input Parameter:
2190: . pc - the preconditioner context

2192:   Output Parameter:
2193: . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

2195:   Level: Intermediate

2197: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2198: .seealso: PCCompositeSetType()
2199: @*/
2200: PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2201: {
2202:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2207:   *type = jac->type;
2208:   return(0);
2209: }

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

2214:    Logically Collective

2216:    Input Parameters:
2217: +  pc   - the preconditioner context
2218: -  flg  - boolean indicating whether to use field splits defined by the DM

2220:    Options Database Key:
2221: .  -pc_fieldsplit_dm_splits

2223:    Level: Intermediate

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

2227: .seealso: PCFieldSplitGetDMSplits()

2229: @*/
2230: PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
2231: {
2232:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2233:   PetscBool      isfs;

2239:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2240:   if (isfs) {
2241:     jac->dm_splits = flg;
2242:   }
2243:   return(0);
2244: }


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

2250:    Logically Collective

2252:    Input Parameter:
2253: .  pc   - the preconditioner context

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

2258:    Level: Intermediate

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

2262: .seealso: PCFieldSplitSetDMSplits()

2264: @*/
2265: PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
2266: {
2267:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2268:   PetscBool      isfs;

2274:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2275:   if (isfs) {
2276:     if(flg) *flg = jac->dm_splits;
2277:   }
2278:   return(0);
2279: }

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

2284:    Logically Collective

2286:    Input Parameter:
2287: .  pc   - the preconditioner context

2289:    Output Parameter:
2290: .  flg  - boolean indicating whether to detect fields or not

2292:    Level: Intermediate

2294: .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint()

2296: @*/
2297: PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC pc,PetscBool *flg)
2298: {
2299:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

2302:   *flg = jac->detect;
2303:   return(0);
2304: }

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

2309:    Logically Collective

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

2314:    Input Parameter:
2315: .  pc   - the preconditioner context

2317:    Output Parameter:
2318: .  flg  - boolean indicating whether to detect fields or not

2320:    Options Database Key:
2321: .  -pc_fieldsplit_detect_saddle_point

2323:    Level: Intermediate

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

2327: @*/
2328: PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc,PetscBool flg)
2329: {
2330:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

2334:   jac->detect = flg;
2335:   if (jac->detect) {
2336:     PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
2337:     PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELF,NULL);
2338:   }
2339:   return(0);
2340: }

2342: /* -------------------------------------------------------------------------------------*/
2343: /*MC
2344:    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2345:                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.

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

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

2353:    Level: intermediate

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

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

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

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

2376: $     For the Schur complement preconditioner if J = ( A00 A01 )
2377: $                                                    ( A10 A11 )
2378: $     the preconditioner using full factorization is
2379: $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2380: $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
2381:      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
2382: $              S = A11 - A10 ksp(A00) A01
2383:      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
2384:      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2385:      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2386:      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.

2388:      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
2389:      diag gives
2390: $              ( inv(A00)     0   )
2391: $              (   0      -ksp(S) )
2392:      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
2393:      can be turned off with PCFieldSplitSetSchurScale() or by command line -pc_fieldsplit_schur_scale 1.0. The lower factorization is the inverse of
2394: $              (  A00   0 )
2395: $              (  A10   S )
2396:      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2397: $              ( A00 A01 )
2398: $              (  0   S  )
2399:      where again the inverses of A00 and S are applied using KSPs.

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

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

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

2411:    Concepts: physics based preconditioners, block preconditioners

2413:    There is a nice discussion of block preconditioners in

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

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

2422: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
2423:            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
2424:           MatSchurComplementSetAinvType(), PCFieldSplitSetSchurScale(),
2425:           PCFieldSplitSetDetectSaddlePoint()
2426: M*/

2428: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
2429: {
2431:   PC_FieldSplit  *jac;

2434:   PetscNewLog(pc,&jac);

2436:   jac->bs                 = -1;
2437:   jac->nsplits            = 0;
2438:   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2439:   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2440:   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2441:   jac->schurscale         = -1.0;
2442:   jac->dm_splits          = PETSC_TRUE;
2443:   jac->detect             = PETSC_FALSE;

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

2447:   pc->ops->apply           = PCApply_FieldSplit;
2448:   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
2449:   pc->ops->setup           = PCSetUp_FieldSplit;
2450:   pc->ops->reset           = PCReset_FieldSplit;
2451:   pc->ops->destroy         = PCDestroy_FieldSplit;
2452:   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
2453:   pc->ops->view            = PCView_FieldSplit;
2454:   pc->ops->applyrichardson = 0;

2456:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2457:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);
2458:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);
2459:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);
2460:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);
2461:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);
2462:   return(0);
2463: }