Actual source code: fieldsplit.c

petsc-master 2015-03-27
Report Typos and Errors
  2: #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
  3: #include <petscdm.h>


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

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

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

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

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

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


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

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

125:  if (isdraw) {
126:     PetscDraw draw;
127:     PetscReal x,y,w,wd;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

494:     jac->issetup = PETSC_TRUE;

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

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

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

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

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

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

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

635:   if (jac->type == PC_COMPOSITE_SCHUR) {
636:     IS          ccis;
637:     PetscInt    rstart,rend;
638:     char        lscname[256];
639:     PetscObject LSC_L;

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

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

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

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

686:       /* extract the A01 and A10 matrices */
687:       ilink = jac->head;
688:       ISComplement(ilink->is_col,rstart,rend,&ccis);
689:       if (jac->offdiag_use_amat) {
690:         MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
691:       } else {
692:         MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
693:       }
694:       ISDestroy(&ccis);
695:       ilink = ilink->next;
696:       ISComplement(ilink->is_col,rstart,rend,&ccis);
697:       if (jac->offdiag_use_amat) {
698:         MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
699:       } else {
700:         MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
701:       }
702:       ISDestroy(&ccis);

704:       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
705:       MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);
706:       MatSetType(jac->schur,MATSCHURCOMPLEMENT);
707:       MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
708:       PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
709:       /* 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? */
710:       MatSetOptionsPrefix(jac->schur,schurmatprefix);
711:       MatSetFromOptions(jac->schur);
712:       MatGetNullSpace(jac->pmat[1], &sp);
713:       if (sp) {
714:         MatSetNullSpace(jac->schur, sp);
715:       }

717:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);
718:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
719:       if (flg) {
720:         DM  dmInner;
721:         KSP kspInner;

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

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

747:       PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);
748:       PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
749:       if (flg) {
750:         DM dmInner;

752:         PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
753:         KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);
754:         KSPSetOptionsPrefix(jac->kspupper, schurprefix);
755:         KSPGetDM(jac->head->ksp, &dmInner);
756:         KSPSetDM(jac->kspupper, dmInner);
757:         KSPSetDMActive(jac->kspupper, PETSC_FALSE);
758:         KSPSetFromOptions(jac->kspupper);
759:         KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);
760:         VecDuplicate(jac->head->x, &jac->head->z);
761:       } else {
762:         jac->kspupper = jac->head->ksp;
763:         PetscObjectReference((PetscObject) jac->head->ksp);
764:       }

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

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

819:   jac->suboptionsset = PETSC_TRUE;
820:   return(0);
821: }

823: #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
824:   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
825:    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
826:    KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
827:    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
828:    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))

832: static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
833: {
834:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
835:   PetscErrorCode    ierr;
836:   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
837:   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;

840:   switch (jac->schurfactorization) {
841:   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
842:     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
843:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
844:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
845:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
846:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
847:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
848:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
849:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
850:     VecScale(ilinkD->y,-1.);
851:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
852:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
853:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
854:     break;
855:   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
856:     /* [A00 0; A10 S], suitable for left preconditioning */
857:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
858:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
859:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
860:     MatMult(jac->C,ilinkA->y,ilinkD->x);
861:     VecScale(ilinkD->x,-1.);
862:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
863:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
864:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
865:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
866:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
867:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
868:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
869:     break;
870:   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
871:     /* [A00 A01; 0 S], suitable for right preconditioning */
872:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
873:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
874:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
875:     MatMult(jac->B,ilinkD->y,ilinkA->x);
876:     VecScale(ilinkA->x,-1.);
877:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
878:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
879:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
880:     KSPSolve(kspA,ilinkA->x,ilinkA->y);
881:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
882:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
883:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
884:     break;
885:   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
886:     /* [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 */
887:     VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
888:     VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
889:     KSPSolve(kspLower,ilinkA->x,ilinkA->y);
890:     MatMult(jac->C,ilinkA->y,ilinkD->x);
891:     VecScale(ilinkD->x,-1.0);
892:     VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
893:     VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);

895:     KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
896:     VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
897:     VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);

899:     if (kspUpper == kspA) {
900:       MatMult(jac->B,ilinkD->y,ilinkA->y);
901:       VecAXPY(ilinkA->x,-1.0,ilinkA->y);
902:       KSPSolve(kspA,ilinkA->x,ilinkA->y);
903:     } else {
904:       KSPSolve(kspA,ilinkA->x,ilinkA->y);
905:       MatMult(jac->B,ilinkD->y,ilinkA->x);
906:       KSPSolve(kspUpper,ilinkA->x,ilinkA->z);
907:       VecAXPY(ilinkA->y,-1.0,ilinkA->z);
908:     }
909:     VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
910:     VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
911:   }
912:   return(0);
913: }

917: static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
918: {
919:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
920:   PetscErrorCode    ierr;
921:   PC_FieldSplitLink ilink = jac->head;
922:   PetscInt          cnt,bs;

925:   if (jac->type == PC_COMPOSITE_ADDITIVE) {
926:     if (jac->defaultsplit) {
927:       VecGetBlockSize(x,&bs);
928:       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);
929:       VecGetBlockSize(y,&bs);
930:       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);
931:       VecStrideGatherAll(x,jac->x,INSERT_VALUES);
932:       while (ilink) {
933:         KSPSolve(ilink->ksp,ilink->x,ilink->y);
934:         ilink = ilink->next;
935:       }
936:       VecStrideScatterAll(jac->y,y,INSERT_VALUES);
937:     } else {
938:       VecSet(y,0.0);
939:       while (ilink) {
940:         FieldSplitSplitSolveAdd(ilink,x,y);
941:         ilink = ilink->next;
942:       }
943:     }
944:   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
945:     VecSet(y,0.0);
946:     /* solve on first block for first block variables */
947:     VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
948:     VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
949:     KSPSolve(ilink->ksp,ilink->x,ilink->y);
950:     VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
951:     VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);

953:     /* compute the residual only onto second block variables using first block variables */
954:     MatMult(jac->Afield[1],ilink->y,ilink->next->x);
955:     ilink = ilink->next;
956:     VecScale(ilink->x,-1.0);
957:     VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
958:     VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);

960:     /* solve on second block variables */
961:     KSPSolve(ilink->ksp,ilink->x,ilink->y);
962:     VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
963:     VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
964:   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
965:     if (!jac->w1) {
966:       VecDuplicate(x,&jac->w1);
967:       VecDuplicate(x,&jac->w2);
968:     }
969:     VecSet(y,0.0);
970:     FieldSplitSplitSolveAdd(ilink,x,y);
971:     cnt  = 1;
972:     while (ilink->next) {
973:       ilink = ilink->next;
974:       /* compute the residual only over the part of the vector needed */
975:       MatMult(jac->Afield[cnt++],y,ilink->x);
976:       VecScale(ilink->x,-1.0);
977:       VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
978:       VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
979:       KSPSolve(ilink->ksp,ilink->x,ilink->y);
980:       VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
981:       VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
982:     }
983:     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
984:       cnt -= 2;
985:       while (ilink->previous) {
986:         ilink = ilink->previous;
987:         /* compute the residual only over the part of the vector needed */
988:         MatMult(jac->Afield[cnt--],y,ilink->x);
989:         VecScale(ilink->x,-1.0);
990:         VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
991:         VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
992:         KSPSolve(ilink->ksp,ilink->x,ilink->y);
993:         VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
994:         VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
995:       }
996:     }
997:   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
998:   return(0);
999: }

1001: #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1002:   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1003:    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1004:    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
1005:    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1006:    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))

1010: static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1011: {
1012:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1013:   PetscErrorCode    ierr;
1014:   PC_FieldSplitLink ilink = jac->head;
1015:   PetscInt          bs;

1018:   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1019:     if (jac->defaultsplit) {
1020:       VecGetBlockSize(x,&bs);
1021:       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);
1022:       VecGetBlockSize(y,&bs);
1023:       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);
1024:       VecStrideGatherAll(x,jac->x,INSERT_VALUES);
1025:       while (ilink) {
1026:         KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);
1027:         ilink = ilink->next;
1028:       }
1029:       VecStrideScatterAll(jac->y,y,INSERT_VALUES);
1030:     } else {
1031:       VecSet(y,0.0);
1032:       while (ilink) {
1033:         FieldSplitSplitSolveAddTranspose(ilink,x,y);
1034:         ilink = ilink->next;
1035:       }
1036:     }
1037:   } else {
1038:     if (!jac->w1) {
1039:       VecDuplicate(x,&jac->w1);
1040:       VecDuplicate(x,&jac->w2);
1041:     }
1042:     VecSet(y,0.0);
1043:     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1044:       FieldSplitSplitSolveAddTranspose(ilink,x,y);
1045:       while (ilink->next) {
1046:         ilink = ilink->next;
1047:         MatMultTranspose(pc->mat,y,jac->w1);
1048:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1049:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1050:       }
1051:       while (ilink->previous) {
1052:         ilink = ilink->previous;
1053:         MatMultTranspose(pc->mat,y,jac->w1);
1054:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1055:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1056:       }
1057:     } else {
1058:       while (ilink->next) {   /* get to last entry in linked list */
1059:         ilink = ilink->next;
1060:       }
1061:       FieldSplitSplitSolveAddTranspose(ilink,x,y);
1062:       while (ilink->previous) {
1063:         ilink = ilink->previous;
1064:         MatMultTranspose(pc->mat,y,jac->w1);
1065:         VecWAXPY(jac->w2,-1.0,jac->w1,x);
1066:         FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1067:       }
1068:     }
1069:   }
1070:   return(0);
1071: }

1075: static PetscErrorCode PCReset_FieldSplit(PC pc)
1076: {
1077:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1078:   PetscErrorCode    ierr;
1079:   PC_FieldSplitLink ilink = jac->head,next;

1082:   while (ilink) {
1083:     KSPReset(ilink->ksp);
1084:     VecDestroy(&ilink->x);
1085:     VecDestroy(&ilink->y);
1086:     VecDestroy(&ilink->z);
1087:     VecScatterDestroy(&ilink->sctx);
1088:     ISDestroy(&ilink->is);
1089:     ISDestroy(&ilink->is_col);
1090:     next  = ilink->next;
1091:     ilink = next;
1092:   }
1093:   PetscFree2(jac->x,jac->y);
1094:   if (jac->mat && jac->mat != jac->pmat) {
1095:     MatDestroyMatrices(jac->nsplits,&jac->mat);
1096:   } else if (jac->mat) {
1097:     jac->mat = NULL;
1098:   }
1099:   if (jac->pmat) {MatDestroyMatrices(jac->nsplits,&jac->pmat);}
1100:   if (jac->Afield) {MatDestroyMatrices(jac->nsplits,&jac->Afield);}
1101:   VecDestroy(&jac->w1);
1102:   VecDestroy(&jac->w2);
1103:   MatDestroy(&jac->schur);
1104:   MatDestroy(&jac->schurp);
1105:   MatDestroy(&jac->schur_user);
1106:   KSPDestroy(&jac->kspschur);
1107:   KSPDestroy(&jac->kspupper);
1108:   MatDestroy(&jac->B);
1109:   MatDestroy(&jac->C);
1110:   jac->reset = PETSC_TRUE;
1111:   return(0);
1112: }

1116: static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1117: {
1118:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1119:   PetscErrorCode    ierr;
1120:   PC_FieldSplitLink ilink = jac->head,next;

1123:   PCReset_FieldSplit(pc);
1124:   while (ilink) {
1125:     KSPDestroy(&ilink->ksp);
1126:     next  = ilink->next;
1127:     PetscFree(ilink->splitname);
1128:     PetscFree(ilink->fields);
1129:     PetscFree(ilink->fields_col);
1130:     PetscFree(ilink);
1131:     ilink = next;
1132:   }
1133:   PetscFree2(jac->x,jac->y);
1134:   PetscFree(pc->data);
1135:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);
1136:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);
1137:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);
1138:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);
1139:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);
1140:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);
1141:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);
1142:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);
1143:   return(0);
1144: }

1148: static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptions *PetscOptionsObject,PC pc)
1149: {
1150:   PetscErrorCode  ierr;
1151:   PetscInt        bs;
1152:   PetscBool       flg,stokes = PETSC_FALSE;
1153:   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
1154:   PCCompositeType ctype;

1157:   PetscOptionsHead(PetscOptionsObject,"FieldSplit options");
1158:   PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);
1159:   PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);
1160:   if (flg) {
1161:     PCFieldSplitSetBlockSize(pc,bs);
1162:   }
1163:   jac->diag_use_amat = pc->useAmat;
1164:   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);
1165:   jac->offdiag_use_amat = pc->useAmat;
1166:   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);
1167:   /* FIXME: No programmatic equivalent to the following. */
1168:   PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);
1169:   if (stokes) {
1170:     PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
1171:     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1172:   }

1174:   PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);
1175:   if (flg) {
1176:     PCFieldSplitSetType(pc,ctype);
1177:   }
1178:   /* Only setup fields once */
1179:   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1180:     /* only allow user to set fields from command line if bs is already known.
1181:        otherwise user can set them in PCFieldSplitSetDefaults() */
1182:     PCFieldSplitSetRuntimeSplits_Private(pc);
1183:     if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
1184:   }
1185:   if (jac->type == PC_COMPOSITE_SCHUR) {
1186:     PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);
1187:     if (flg) {PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");}
1188:     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);
1189:     PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);
1190:   }
1191:   PetscOptionsTail();
1192:   return(0);
1193: }

1195: /*------------------------------------------------------------------------------------*/

1199: static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1200: {
1201:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1202:   PetscErrorCode    ierr;
1203:   PC_FieldSplitLink ilink,next = jac->head;
1204:   char              prefix[128];
1205:   PetscInt          i;

1208:   if (jac->splitdefined) {
1209:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1210:     return(0);
1211:   }
1212:   for (i=0; i<n; i++) {
1213:     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1214:     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1215:   }
1216:   PetscNew(&ilink);
1217:   if (splitname) {
1218:     PetscStrallocpy(splitname,&ilink->splitname);
1219:   } else {
1220:     PetscMalloc1(3,&ilink->splitname);
1221:     PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);
1222:   }
1223:   PetscMalloc1(n,&ilink->fields);
1224:   PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));
1225:   PetscMalloc1(n,&ilink->fields_col);
1226:   PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));

1228:   ilink->nfields = n;
1229:   ilink->next    = NULL;
1230:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1231:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1232:   KSPSetType(ilink->ksp,KSPPREONLY);
1233:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

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

1238:   if (!next) {
1239:     jac->head       = ilink;
1240:     ilink->previous = NULL;
1241:   } else {
1242:     while (next->next) {
1243:       next = next->next;
1244:     }
1245:     next->next      = ilink;
1246:     ilink->previous = next;
1247:   }
1248:   jac->nsplits++;
1249:   return(0);
1250: }

1254: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1255: {
1256:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1260:   PetscMalloc1(jac->nsplits,subksp);
1261:   MatSchurComplementGetKSP(jac->schur,*subksp);

1263:   (*subksp)[1] = jac->kspschur;
1264:   if (n) *n = jac->nsplits;
1265:   return(0);
1266: }

1270: static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1271: {
1272:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1273:   PetscErrorCode    ierr;
1274:   PetscInt          cnt   = 0;
1275:   PC_FieldSplitLink ilink = jac->head;

1278:   PetscMalloc1(jac->nsplits,subksp);
1279:   while (ilink) {
1280:     (*subksp)[cnt++] = ilink->ksp;
1281:     ilink            = ilink->next;
1282:   }
1283:   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);
1284:   if (n) *n = jac->nsplits;
1285:   return(0);
1286: }

1290: static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1291: {
1292:   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1293:   PetscErrorCode    ierr;
1294:   PC_FieldSplitLink ilink, next = jac->head;
1295:   char              prefix[128];

1298:   if (jac->splitdefined) {
1299:     PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1300:     return(0);
1301:   }
1302:   PetscNew(&ilink);
1303:   if (splitname) {
1304:     PetscStrallocpy(splitname,&ilink->splitname);
1305:   } else {
1306:     PetscMalloc1(8,&ilink->splitname);
1307:     PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);
1308:   }
1309:   PetscObjectReference((PetscObject)is);
1310:   ISDestroy(&ilink->is);
1311:   ilink->is     = is;
1312:   PetscObjectReference((PetscObject)is);
1313:   ISDestroy(&ilink->is_col);
1314:   ilink->is_col = is;
1315:   ilink->next   = NULL;
1316:   KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1317:   PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1318:   KSPSetType(ilink->ksp,KSPPREONLY);
1319:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);

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

1324:   if (!next) {
1325:     jac->head       = ilink;
1326:     ilink->previous = NULL;
1327:   } else {
1328:     while (next->next) {
1329:       next = next->next;
1330:     }
1331:     next->next      = ilink;
1332:     ilink->previous = next;
1333:   }
1334:   jac->nsplits++;
1335:   return(0);
1336: }

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

1343:     Logically Collective on PC

1345:     Input Parameters:
1346: +   pc  - the preconditioner context
1347: .   splitname - name of this split, if NULL the number of the split is used
1348: .   n - the number of fields in this split
1349: -   fields - the fields in this split

1351:     Level: intermediate

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

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

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

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

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

1369: @*/
1370: PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1371: {

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

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

1388:     Logically Collective on PC

1390:     Input Parameters:
1391: +   pc  - the preconditioner object
1392: -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from

1394:     Options Database:
1395: .     -pc_fieldsplit_diag_use_amat

1397:     Level: intermediate

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

1401: @*/
1402: PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1403: {
1404:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1405:   PetscBool      isfs;

1410:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1411:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1412:   jac->diag_use_amat = flg;
1413:   return(0);
1414: }

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

1421:     Logically Collective on PC

1423:     Input Parameters:
1424: .   pc  - the preconditioner object

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


1430:     Level: intermediate

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

1434: @*/
1435: PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1436: {
1437:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1438:   PetscBool      isfs;

1444:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1445:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1446:   *flg = jac->diag_use_amat;
1447:   return(0);
1448: }

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

1455:     Logically Collective on PC

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

1461:     Options Database:
1462: .     -pc_fieldsplit_off_diag_use_amat

1464:     Level: intermediate

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

1468: @*/
1469: PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1470: {
1471:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1472:   PetscBool      isfs;

1477:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1478:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1479:   jac->offdiag_use_amat = flg;
1480:   return(0);
1481: }

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

1488:     Logically Collective on PC

1490:     Input Parameters:
1491: .   pc  - the preconditioner object

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


1497:     Level: intermediate

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

1501: @*/
1502: PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1503: {
1504:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1505:   PetscBool      isfs;

1511:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1512:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1513:   *flg = jac->offdiag_use_amat;
1514:   return(0);
1515: }



1521: /*@C
1522:     PCFieldSplitSetIS - Sets the exact elements for field

1524:     Logically Collective on PC

1526:     Input Parameters:
1527: +   pc  - the preconditioner context
1528: .   splitname - name of this split, if NULL the number of the split is used
1529: -   is - the index set that defines the vector elements in this field


1532:     Notes:
1533:     Use PCFieldSplitSetFields(), for fields defined by strided types.

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

1538:     Level: intermediate

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

1542: @*/
1543: PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1544: {

1551:   PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));
1552:   return(0);
1553: }

1557: /*@
1558:     PCFieldSplitGetIS - Retrieves the elements for a field as an IS

1560:     Logically Collective on PC

1562:     Input Parameters:
1563: +   pc  - the preconditioner context
1564: -   splitname - name of this split

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

1569:     Level: intermediate

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

1573: @*/
1574: PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
1575: {

1582:   {
1583:     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
1584:     PC_FieldSplitLink ilink = jac->head;
1585:     PetscBool         found;

1587:     *is = NULL;
1588:     while (ilink) {
1589:       PetscStrcmp(ilink->splitname, splitname, &found);
1590:       if (found) {
1591:         *is = ilink->is;
1592:         break;
1593:       }
1594:       ilink = ilink->next;
1595:     }
1596:   }
1597:   return(0);
1598: }

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

1606:     Logically Collective on PC

1608:     Input Parameters:
1609: +   pc  - the preconditioner context
1610: -   bs - the block size

1612:     Level: intermediate

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

1616: @*/
1617: PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1618: {

1624:   PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));
1625:   return(0);
1626: }

1630: /*@C
1631:    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits

1633:    Collective on KSP

1635:    Input Parameter:
1636: .  pc - the preconditioner context

1638:    Output Parameters:
1639: +  n - the number of splits
1640: -  pc - the array of KSP contexts

1642:    Note:
1643:    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1644:    (not the KSP just the array that contains them).

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

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


1653:    Level: advanced

1655: .seealso: PCFIELDSPLIT
1656: @*/
1657: PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1658: {

1664:   PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
1665:   return(0);
1666: }

1670: /*@
1671:     PCFieldSplitSetSchurPre -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1672:       A11 matrix. Otherwise no preconditioner is used.

1674:     Collective on PC

1676:     Input Parameters:
1677: +   pc      - the preconditioner context
1678: .   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
1679: -   userpre - matrix to use for preconditioning, or NULL

1681:     Options Database:
1682: .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> default is a11

1684:     Notes:
1685: $    If ptype is
1686: $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1687: $             to this function).
1688: $        a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the preconditioner
1689: $             matrix associated with the Schur complement (i.e. A11), not he Schur complement matrix
1690: $        full then the preconditioner uses the exact Schur complement (this is expensive)
1691: $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1692: $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1693: $             preconditioner
1694: $        selfp then the preconditioning matrix is an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
1695: $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00.

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

1701:     Level: intermediate

1703: .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC

1705: @*/
1706: PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1707: {

1712:   PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));
1713:   return(0);
1714: }
1715: PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */

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

1723:     Logically Collective on PC

1725:     Input Parameters:
1726: .   pc      - the preconditioner context

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

1732:     Level: intermediate

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

1736: @*/
1737: PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
1738: {

1743:   PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));
1744:   return(0);
1745: }

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

1752:     Not collective

1754:     Input Parameter:
1755: .   pc      - the preconditioner context

1757:     Output Parameter:
1758: .   S       - the Schur complement matrix

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

1763:     Level: advanced

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

1767: @*/
1768: PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
1769: {
1771:   const char*    t;
1772:   PetscBool      isfs;
1773:   PC_FieldSplit  *jac;

1777:   PetscObjectGetType((PetscObject)pc,&t);
1778:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
1779:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1780:   jac = (PC_FieldSplit*)pc->data;
1781:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1782:   if (S) *S = jac->schur;
1783:   return(0);
1784: }

1788: /*@
1789:     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC

1791:     Not collective

1793:     Input Parameters:
1794: +   pc      - the preconditioner context
1795: .   S       - the Schur complement matrix

1797:     Level: advanced

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

1801: @*/
1802: PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
1803: {
1805:   const char*    t;
1806:   PetscBool      isfs;
1807:   PC_FieldSplit  *jac;

1811:   PetscObjectGetType((PetscObject)pc,&t);
1812:   PetscStrcmp(t,PCFIELDSPLIT,&isfs);
1813:   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1814:   jac = (PC_FieldSplit*)pc->data;
1815:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1816:   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
1817:   return(0);
1818: }


1823: static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1824: {
1825:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1829:   jac->schurpre = ptype;
1830:   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
1831:     MatDestroy(&jac->schur_user);
1832:     jac->schur_user = pre;
1833:     PetscObjectReference((PetscObject)jac->schur_user);
1834:   }
1835:   return(0);
1836: }

1840: static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
1841: {
1842:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1845:   *ptype = jac->schurpre;
1846:   *pre   = jac->schur_user;
1847:   return(0);
1848: }

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

1855:     Collective on PC

1857:     Input Parameters:
1858: +   pc  - the preconditioner context
1859: -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default

1861:     Options Database:
1862: .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full


1865:     Level: intermediate

1867:     Notes:
1868:     The FULL factorization is

1870: $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1871: $   (C   D)    (C*Ainv  1) (0   S) (0     1  )

1873:     where S = D - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D,
1874:     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).

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

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

1884:     References:
1885:     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.

1887:     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.

1889: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1890: @*/
1891: PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1892: {

1897:   PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));
1898:   return(0);
1899: }

1903: static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1904: {
1905:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

1908:   jac->schurfactorization = ftype;
1909:   return(0);
1910: }

1914: /*@C
1915:    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement

1917:    Collective on KSP

1919:    Input Parameter:
1920: .  pc - the preconditioner context

1922:    Output Parameters:
1923: +  A00 - the (0,0) block
1924: .  A01 - the (0,1) block
1925: .  A10 - the (1,0) block
1926: -  A11 - the (1,1) block

1928:    Level: advanced

1930: .seealso: PCFIELDSPLIT
1931: @*/
1932: PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
1933: {
1934:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

1938:   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1939:   if (A00) *A00 = jac->pmat[0];
1940:   if (A01) *A01 = jac->B;
1941:   if (A10) *A10 = jac->C;
1942:   if (A11) *A11 = jac->pmat[1];
1943:   return(0);
1944: }

1948: static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
1949: {
1950:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;

1954:   jac->type = type;
1955:   if (type == PC_COMPOSITE_SCHUR) {
1956:     pc->ops->apply = PCApply_FieldSplit_Schur;
1957:     pc->ops->view  = PCView_FieldSplit_Schur;

1959:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);
1960:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);
1961:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);
1962:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);

1964:   } else {
1965:     pc->ops->apply = PCApply_FieldSplit;
1966:     pc->ops->view  = PCView_FieldSplit;

1968:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
1969:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);
1970:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);
1971:     PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);
1972:   }
1973:   return(0);
1974: }

1978: static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
1979: {
1980:   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;

1983:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1984:   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);
1985:   jac->bs = bs;
1986:   return(0);
1987: }

1991: /*@
1992:    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.

1994:    Collective on PC

1996:    Input Parameter:
1997: .  pc - the preconditioner context
1998: .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

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

2003:    Level: Intermediate

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

2007: .seealso: PCCompositeSetType()

2009: @*/
2010: PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
2011: {

2016:   PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));
2017:   return(0);
2018: }

2022: /*@
2023:   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.

2025:   Not collective

2027:   Input Parameter:
2028: . pc - the preconditioner context

2030:   Output Parameter:
2031: . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

2033:   Level: Intermediate

2035: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2036: .seealso: PCCompositeSetType()
2037: @*/
2038: PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2039: {
2040:   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;

2045:   *type = jac->type;
2046:   return(0);
2047: }

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

2054:    Logically Collective

2056:    Input Parameters:
2057: +  pc   - the preconditioner context
2058: -  flg  - boolean indicating whether to use field splits defined by the DM

2060:    Options Database Key:
2061: .  -pc_fieldsplit_dm_splits

2063:    Level: Intermediate

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

2067: .seealso: PCFieldSplitGetDMSplits()

2069: @*/
2070: PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
2071: {
2072:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2073:   PetscBool      isfs;

2079:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2080:   if (isfs) {
2081:     jac->dm_splits = flg;
2082:   }
2083:   return(0);
2084: }


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

2092:    Logically Collective

2094:    Input Parameter:
2095: .  pc   - the preconditioner context

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

2100:    Level: Intermediate

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

2104: .seealso: PCFieldSplitSetDMSplits()

2106: @*/
2107: PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
2108: {
2109:   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2110:   PetscBool      isfs;

2116:   PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2117:   if (isfs) {
2118:     if(flg) *flg = jac->dm_splits;
2119:   }
2120:   return(0);
2121: }

2123: /* -------------------------------------------------------------------------------------*/
2124: /*MC
2125:    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2126:                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.

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

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

2134:    Level: intermediate

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

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

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

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

2157: $     For the Schur complement preconditioner if J = ( A00 A01 )
2158: $                                                    ( A10 A11 )
2159: $     the preconditioner using full factorization is
2160: $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2161: $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
2162:      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
2163: $              S = A11 - A10 ksp(A00) A01
2164:      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
2165:      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2166:      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2167:      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() to turn on or off this
2168:      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc.
2169:      When option -fieldsplit_schur_precondition selfp is given, an approximation to S is assembled --
2170:      Sp = A11 - A10 inv(diag(A00)) A01, which has type AIJ and can be used with a variety of preconditioners
2171:      (e.g., -fieldsplit_1_pc_type asm).
2172:      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
2173:      diag gives
2174: $              ( inv(A00)     0   )
2175: $              (   0      -ksp(S) )
2176:      note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of
2177: $              (  A00   0 )
2178: $              (  A10   S )
2179:      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2180: $              ( A00 A01 )
2181: $              (  0   S  )
2182:      where again the inverses of A00 and S are applied using KSPs.

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

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

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

2194:    Concepts: physics based preconditioners, block preconditioners

2196:    There is a nice discussion of block preconditioners in

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

2202: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
2203:            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre()
2204: M*/

2208: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
2209: {
2211:   PC_FieldSplit  *jac;

2214:   PetscNewLog(pc,&jac);

2216:   jac->bs                 = -1;
2217:   jac->nsplits            = 0;
2218:   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2219:   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2220:   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2221:   jac->dm_splits          = PETSC_TRUE;

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

2225:   pc->ops->apply           = PCApply_FieldSplit;
2226:   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
2227:   pc->ops->setup           = PCSetUp_FieldSplit;
2228:   pc->ops->reset           = PCReset_FieldSplit;
2229:   pc->ops->destroy         = PCDestroy_FieldSplit;
2230:   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
2231:   pc->ops->view            = PCView_FieldSplit;
2232:   pc->ops->applyrichardson = 0;

2234:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2235:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);
2236:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);
2237:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);
2238:   PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);
2239:   return(0);
2240: }