Actual source code: fieldsplit.c

petsc-master 2015-03-03
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:     PetscDrawBoxedString(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:     PetscDrawBoxedString(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:         } else {
521:           ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);
522:           ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);
523:        }
524:       }
525:       ISSorted(ilink->is,&sorted);
526:       if (ilink->is_col) { ISSorted(ilink->is_col,&sorted_col); }
527:       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
528:       ilink = ilink->next;
529:     }
530:   }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

817:   jac->suboptionsset = PETSC_TRUE;
818:   return(0);
819: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1193: /*------------------------------------------------------------------------------------*/

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

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

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

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

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

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

1258:   PetscMalloc1(jac->nsplits,subksp);
1259:   MatSchurComplementGetKSP(jac->schur,*subksp);

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

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

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

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

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

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

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

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

1341:     Logically Collective on PC

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

1349:     Level: intermediate

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

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

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

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

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

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

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

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

1386:     Logically Collective on PC

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

1392:     Options Database:
1393: .     -pc_fieldsplit_diag_use_amat

1395:     Level: intermediate

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

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

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

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

1419:     Logically Collective on PC

1421:     Input Parameters:
1422: .   pc  - the preconditioner object

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


1428:     Level: intermediate

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

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

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

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

1453:     Logically Collective on PC

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

1459:     Options Database:
1460: .     -pc_fieldsplit_off_diag_use_amat

1462:     Level: intermediate

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

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

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

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

1486:     Logically Collective on PC

1488:     Input Parameters:
1489: .   pc  - the preconditioner object

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


1495:     Level: intermediate

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

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

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



1519: /*@C
1520:     PCFieldSplitSetIS - Sets the exact elements for field

1522:     Logically Collective on PC

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


1530:     Notes:
1531:     Use PCFieldSplitSetFields(), for fields defined by strided types.

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

1536:     Level: intermediate

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

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

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

1555: /*@
1556:     PCFieldSplitGetIS - Retrieves the elements for a field as an IS

1558:     Logically Collective on PC

1560:     Input Parameters:
1561: +   pc  - the preconditioner context
1562: -   splitname - name of this split

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

1567:     Level: intermediate

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

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

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

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

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

1604:     Logically Collective on PC

1606:     Input Parameters:
1607: +   pc  - the preconditioner context
1608: -   bs - the block size

1610:     Level: intermediate

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

1614: @*/
1615: PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1616: {

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

1628: /*@C
1629:    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits

1631:    Collective on KSP

1633:    Input Parameter:
1634: .  pc - the preconditioner context

1636:    Output Parameters:
1637: +  n - the number of splits
1638: -  pc - the array of KSP contexts

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

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

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


1651:    Level: advanced

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

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

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

1672:     Collective on PC

1674:     Input Parameters:
1675: +   pc      - the preconditioner context
1676: .   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
1677: -   userpre - matrix to use for preconditioning, or NULL

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

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

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

1699:     Level: intermediate

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

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

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

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

1721:     Logically Collective on PC

1723:     Input Parameters:
1724: .   pc      - the preconditioner context

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

1730:     Level: intermediate

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

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

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

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

1750:     Not collective

1752:     Input Parameter:
1753: .   pc      - the preconditioner context

1755:     Output Parameter:
1756: .   S       - the Schur complement matrix

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

1761:     Level: advanced

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

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

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

1786: /*@
1787:     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC

1789:     Not collective

1791:     Input Parameters:
1792: +   pc      - the preconditioner context
1793: .   S       - the Schur complement matrix

1795:     Level: advanced

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

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

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


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

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

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

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

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

1853:     Collective on PC

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

1859:     Options Database:
1860: .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full


1863:     Level: intermediate

1865:     Notes:
1866:     The FULL factorization is

1868: $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1869: $   (C   D)    (C*Ainv  1) (0   S) (0     1  )

1871:     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,
1872:     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).

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

1879:     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
1880:     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).

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

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

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

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

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

1906:   jac->schurfactorization = ftype;
1907:   return(0);
1908: }

1912: /*@C
1913:    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement

1915:    Collective on KSP

1917:    Input Parameter:
1918: .  pc - the preconditioner context

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

1926:    Level: advanced

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

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

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

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

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

1962:   } else {
1963:     pc->ops->apply = PCApply_FieldSplit;
1964:     pc->ops->view  = PCView_FieldSplit;

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

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

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

1989: /*@
1990:    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.

1992:    Collective on PC

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

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

2001:    Level: Intermediate

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

2005: .seealso: PCCompositeSetType()

2007: @*/
2008: PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
2009: {

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

2020: /*@
2021:   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.

2023:   Not collective

2025:   Input Parameter:
2026: . pc - the preconditioner context

2028:   Output Parameter:
2029: . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR

2031:   Level: Intermediate

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

2043:   *type = jac->type;
2044:   return(0);
2045: }

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

2052:    Logically Collective

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

2058:    Options Database Key:
2059: .  -pc_fieldsplit_dm_splits

2061:    Level: Intermediate

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

2065: .seealso: PCFieldSplitGetDMSplits()

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

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


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

2090:    Logically Collective

2092:    Input Parameter:
2093: .  pc   - the preconditioner context

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

2098:    Level: Intermediate

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

2102: .seealso: PCFieldSplitSetDMSplits()

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

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

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

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

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

2132:    Level: intermediate

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

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

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

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

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

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

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

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

2192:    Concepts: physics based preconditioners, block preconditioners

2194:    There is a nice discussion of block preconditioners in

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

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

2206: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
2207: {
2209:   PC_FieldSplit  *jac;

2212:   PetscNewLog(pc,&jac);

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

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

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

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