xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision bdf89e91df64bd49f4414ba801fb07e6338da027)
1 
2 #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
3 #include <petscdm.h>
4 
5 /*
6   There is a nice discussion of block preconditioners in
7 
8 [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier–Stokes equations
9        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
10        http://chess.cs.umd.edu/~elman/papers/tax.pdf
11 */
12 
13 const char *const PCFieldSplitSchurPreTypes[] = {"SELF","A11","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
14 const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
15 
16 typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
17 struct _PC_FieldSplitLink {
18   KSP               ksp;
19   Vec               x,y,z;
20   char              *splitname;
21   PetscInt          nfields;
22   PetscInt          *fields,*fields_col;
23   VecScatter        sctx;
24   IS                is,is_col;
25   PC_FieldSplitLink next,previous;
26 };
27 
28 typedef struct {
29   PCCompositeType type;
30   PetscBool       defaultsplit;                    /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
31   PetscBool       splitdefined;                    /* Flag is set after the splits have been defined, to prevent more splits from being added */
32   PetscInt        bs;                              /* Block size for IS and Mat structures */
33   PetscInt        nsplits;                         /* Number of field divisions defined */
34   Vec             *x,*y,w1,w2;
35   Mat             *mat;                            /* The diagonal block for each split */
36   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
37   Mat             *Afield;                         /* The rows of the matrix associated with each split */
38   PetscBool       issetup;
39 
40   /* Only used when Schur complement preconditioning is used */
41   Mat                       B;                     /* The (0,1) block */
42   Mat                       C;                     /* The (1,0) block */
43   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
44   Mat                       schur_user;            /* User-provided preconditioning matrix for the Schur complement */
45   PCFieldSplitSchurPreType  schurpre;              /* Determines which preconditioning matrix is used for the Schur complement */
46   PCFieldSplitSchurFactType schurfactorization;
47   KSP                       kspschur;              /* The solver for S */
48   KSP                       kspupper;              /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
49   PC_FieldSplitLink         head;
50   PetscBool                 reset;                  /* indicates PCReset() has been last called on this object, hack */
51   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
52   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
53 } PC_FieldSplit;
54 
55 /*
56     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
57    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
58    PC you could change this.
59 */
60 
61 /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
62 * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
63 static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
64 {
65   switch (jac->schurpre) {
66   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
67   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
68   case PC_FIELDSPLIT_SCHUR_PRE_USER:   /* Use a user-provided matrix if it is given, otherwise diagonal block */
69   default:
70     return jac->schur_user ? jac->schur_user : jac->pmat[1];
71   }
72 }
73 
74 
75 #include <petscdraw.h>
76 #undef __FUNCT__
77 #define __FUNCT__ "PCView_FieldSplit"
78 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
79 {
80   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
81   PetscErrorCode    ierr;
82   PetscBool         iascii,isdraw;
83   PetscInt          i,j;
84   PC_FieldSplitLink ilink = jac->head;
85 
86   PetscFunctionBegin;
87   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
88   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
89   if (iascii) {
90     if (jac->bs > 0) {
91       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
92     } else {
93       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
94     }
95     if (pc->useAmat) {
96       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
97     }
98     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
99     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
100     for (i=0; i<jac->nsplits; i++) {
101       if (ilink->fields) {
102         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
103         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
104         for (j=0; j<ilink->nfields; j++) {
105           if (j > 0) {
106             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
107           }
108           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
109         }
110         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
111         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
112       } else {
113         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
114       }
115       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
116       ilink = ilink->next;
117     }
118     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
119   }
120 
121  if (isdraw) {
122     PetscDraw draw;
123     PetscReal x,y,w,wd;
124 
125     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
126     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
127     w    = 2*PetscMin(1.0 - x,x);
128     wd   = w/(jac->nsplits + 1);
129     x    = x - wd*(jac->nsplits-1)/2.0;
130     for (i=0; i<jac->nsplits; i++) {
131       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
132       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
133       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
134       x    += wd;
135       ilink = ilink->next;
136     }
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "PCView_FieldSplit_Schur"
143 static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
144 {
145   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
146   PetscErrorCode    ierr;
147   PetscBool         iascii,isdraw;
148   PetscInt          i,j;
149   PC_FieldSplitLink ilink = jac->head;
150 
151   PetscFunctionBegin;
152   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
153   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
154   if (iascii) {
155     if (jac->bs > 0) {
156       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
157     } else {
158       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
159     }
160     if (pc->useAmat) {
161       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
162     }
163     switch (jac->schurpre) {
164     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
165       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break;
166     case PC_FIELDSPLIT_SCHUR_PRE_A11:
167       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break;
168     case PC_FIELDSPLIT_SCHUR_PRE_USER:
169       if (jac->schur_user) {
170         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
171       } else {
172         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
173       }
174       break;
175     default:
176       SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
177     }
178     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
179     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
180     for (i=0; i<jac->nsplits; i++) {
181       if (ilink->fields) {
182         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
183         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
184         for (j=0; j<ilink->nfields; j++) {
185           if (j > 0) {
186             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
187           }
188           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
189         }
190         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
191         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
192       } else {
193         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
194       }
195       ilink = ilink->next;
196     }
197     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
198     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
199     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
200     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
201     if (jac->kspupper != jac->head->ksp) {
202       ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr);
203       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
204       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
205       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
206     }
207     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
208     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
209     if (jac->kspschur) {
210       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
211     } else {
212       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
213     }
214     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
215     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
216   } else if (isdraw) {
217     PetscDraw draw;
218     PetscReal x,y,w,wd,h;
219     PetscInt  cnt = 2;
220     char      str[32];
221 
222     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
223     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
224     if (jac->kspupper != jac->head->ksp) cnt++;
225     w  = 2*PetscMin(1.0 - x,x);
226     wd = w/(cnt + 1);
227 
228     ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
229     ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
230     y   -= h;
231     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
232       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr);
233     } else {
234       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr);
235     }
236     ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
237     y   -= h;
238     x    = x - wd*(cnt-1)/2.0;
239 
240     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
241     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
242     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
243     if (jac->kspupper != jac->head->ksp) {
244       x   += wd;
245       ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
246       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
247       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
248     }
249     x   += wd;
250     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
251     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
252     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
253   }
254   PetscFunctionReturn(0);
255 }
256 
257 #undef __FUNCT__
258 #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
259 /* Precondition: jac->bs is set to a meaningful value */
260 static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
261 {
262   PetscErrorCode ierr;
263   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
264   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
265   PetscBool      flg,flg_col;
266   char           optionname[128],splitname[8],optionname_col[128];
267 
268   PetscFunctionBegin;
269   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
270   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
271   for (i=0,flg=PETSC_TRUE;; i++) {
272     ierr        = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
273     ierr        = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
274     ierr        = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
275     nfields     = jac->bs;
276     nfields_col = jac->bs;
277     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
278     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
279     if (!flg) break;
280     else if (flg && !flg_col) {
281       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
282       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
283     } else {
284       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
285       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
286       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
287     }
288   }
289   if (i > 0) {
290     /* Makes command-line setting of splits take precedence over setting them in code.
291        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
292        create new splits, which would probably not be what the user wanted. */
293     jac->splitdefined = PETSC_TRUE;
294   }
295   ierr = PetscFree(ifields);CHKERRQ(ierr);
296   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "PCFieldSplitSetDefaults"
302 static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
303 {
304   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
305   PetscErrorCode    ierr;
306   PC_FieldSplitLink ilink              = jac->head;
307   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE;
308   PetscInt          i;
309 
310   PetscFunctionBegin;
311   /*
312    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
313    Should probably be rewritten.
314    */
315   if (!ilink || jac->reset) {
316     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
317     if (pc->dm && jac->dm_splits && !stokes) {
318       PetscInt  numFields, f, i, j;
319       char      **fieldNames;
320       IS        *fields;
321       DM        *dms;
322       DM        subdm[128];
323       PetscBool flg;
324 
325       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
326       /* Allow the user to prescribe the splits */
327       for (i = 0, flg = PETSC_TRUE;; i++) {
328         PetscInt ifields[128];
329         IS       compField;
330         char     optionname[128], splitname[8];
331         PetscInt nfields = numFields;
332 
333         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
334         ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
335         if (!flg) break;
336         if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
337         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
338         if (nfields == 1) {
339           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
340           /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr);
341              ierr = ISView(compField, NULL);CHKERRQ(ierr); */
342         } else {
343           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
344           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
345           /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", splitname);CHKERRQ(ierr);
346              ierr = ISView(compField, NULL);CHKERRQ(ierr); */
347         }
348         ierr = ISDestroy(&compField);CHKERRQ(ierr);
349         for (j = 0; j < nfields; ++j) {
350           f    = ifields[j];
351           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
352           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
353         }
354       }
355       if (i == 0) {
356         for (f = 0; f < numFields; ++f) {
357           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
358           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
359           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
360         }
361       } else {
362         for (j=0; j<numFields; j++) {
363           ierr = DMDestroy(dms+j);CHKERRQ(ierr);
364         }
365         ierr = PetscFree(dms);CHKERRQ(ierr);
366         ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr);
367         for (j = 0; j < i; ++j) dms[j] = subdm[j];
368       }
369       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
370       ierr = PetscFree(fields);CHKERRQ(ierr);
371       if (dms) {
372         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
373         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
374           const char *prefix;
375           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr);
376           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr);
377           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
378           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
379           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr);
380           ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
381         }
382         ierr = PetscFree(dms);CHKERRQ(ierr);
383       }
384     } else {
385       if (jac->bs <= 0) {
386         if (pc->pmat) {
387           ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
388         } else jac->bs = 1;
389       }
390 
391       if (stokes) {
392         IS       zerodiags,rest;
393         PetscInt nmin,nmax;
394 
395         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
396         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
397         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
398         if (jac->reset) {
399           jac->head->is       = rest;
400           jac->head->next->is = zerodiags;
401         } else {
402           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
403           ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
404         }
405         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
406         ierr = ISDestroy(&rest);CHKERRQ(ierr);
407       } else {
408         if (jac->reset) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
409         ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
410         if (!fieldsplit_default) {
411           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
412            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
413           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
414           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
415         }
416         if (fieldsplit_default || !jac->splitdefined) {
417           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
418           for (i=0; i<jac->bs; i++) {
419             char splitname[8];
420             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
421             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
422           }
423           jac->defaultsplit = PETSC_TRUE;
424         }
425       }
426     }
427   } else if (jac->nsplits == 1) {
428     if (ilink->is) {
429       IS       is2;
430       PetscInt nmin,nmax;
431 
432       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
433       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
434       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
435       ierr = ISDestroy(&is2);CHKERRQ(ierr);
436     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
437   }
438 
439 
440   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
441   PetscFunctionReturn(0);
442 }
443 
444 PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg);
445 
446 #undef __FUNCT__
447 #define __FUNCT__ "PCSetUp_FieldSplit"
448 static PetscErrorCode PCSetUp_FieldSplit(PC pc)
449 {
450   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
451   PetscErrorCode    ierr;
452   PC_FieldSplitLink ilink;
453   PetscInt          i,nsplit;
454   MatStructure      flag = pc->flag;
455   PetscBool         sorted, sorted_col;
456 
457   PetscFunctionBegin;
458   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
459   nsplit = jac->nsplits;
460   ilink  = jac->head;
461 
462   /* get the matrices for each split */
463   if (!jac->issetup) {
464     PetscInt rstart,rend,nslots,bs;
465 
466     jac->issetup = PETSC_TRUE;
467 
468     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
469     if (jac->defaultsplit || !ilink->is) {
470       if (jac->bs <= 0) jac->bs = nsplit;
471     }
472     bs     = jac->bs;
473     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
474     nslots = (rend - rstart)/bs;
475     for (i=0; i<nsplit; i++) {
476       if (jac->defaultsplit) {
477         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
478         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
479       } else if (!ilink->is) {
480         if (ilink->nfields > 1) {
481           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
482           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
483           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
484           for (j=0; j<nslots; j++) {
485             for (k=0; k<nfields; k++) {
486               ii[nfields*j + k] = rstart + bs*j + fields[k];
487               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
488             }
489           }
490           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
491           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
492         } else {
493           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
494           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
495        }
496       }
497       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
498       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
499       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
500       ilink = ilink->next;
501     }
502   }
503 
504   ilink = jac->head;
505   if (!jac->pmat) {
506     Vec xtmp;
507 
508     ierr = MatGetVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
509     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
510     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
511     for (i=0; i<nsplit; i++) {
512       MatNullSpace sp;
513 
514       /* Check for preconditioning matrix attached to IS */
515       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
516       if (jac->pmat[i]) {
517         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
518         if (jac->type == PC_COMPOSITE_SCHUR) {
519           jac->schur_user = jac->pmat[i];
520 
521           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
522         }
523       } else {
524         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
525       }
526       /* create work vectors for each split */
527       ierr     = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
528       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
529       /* compute scatter contexts needed by multiplicative versions and non-default splits */
530       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
531       /* Check for null space attached to IS */
532       ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
533       if (sp) {
534         ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
535       }
536       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
537       if (sp) {
538         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
539       }
540       ilink = ilink->next;
541     }
542     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
543   } else {
544     for (i=0; i<nsplit; i++) {
545       Mat pmat;
546 
547       /* Check for preconditioning matrix attached to IS */
548       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
549       if (!pmat) {
550         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
551       }
552       ilink = ilink->next;
553     }
554   }
555   if (pc->useAmat) {
556     ilink = jac->head;
557     if (!jac->mat) {
558       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
559       for (i=0; i<nsplit; i++) {
560         ierr  = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
561         ilink = ilink->next;
562       }
563     } else {
564       for (i=0; i<nsplit; i++) {
565         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
566         ilink = ilink->next;
567       }
568     }
569   } else {
570     jac->mat = jac->pmat;
571   }
572 
573   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
574     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
575     ilink = jac->head;
576     if (!jac->Afield) {
577       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
578       for (i=0; i<nsplit; i++) {
579         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
580         ilink = ilink->next;
581       }
582     } else {
583       for (i=0; i<nsplit; i++) {
584         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
585         ilink = ilink->next;
586       }
587     }
588   }
589 
590   if (jac->type == PC_COMPOSITE_SCHUR) {
591     IS          ccis;
592     PetscInt    rstart,rend;
593     char        lscname[256];
594     PetscObject LSC_L;
595 
596     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
597 
598     /* When extracting off-diagonal submatrices, we take complements from this range */
599     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
600 
601     /* need to handle case when one is resetting up the preconditioner */
602     if (jac->schur) {
603       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
604 
605       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
606       ilink = jac->head;
607       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
608       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
609       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
610       ilink = ilink->next;
611       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
612       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
613       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
614       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr);
615       if (kspA != kspInner) {
616         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr);
617       }
618       if (kspUpper != kspA) {
619         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr);
620       }
621       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
622     } else {
623       KSP          ksp;
624       const char   *Dprefix;
625       char         schurprefix[256];
626       char         schurtestoption[256];
627       MatNullSpace sp;
628       PetscBool    flg;
629 
630       /* extract the A01 and A10 matrices */
631       ilink = jac->head;
632       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
633       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
634       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
635       ilink = ilink->next;
636       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
637       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
638       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
639 
640       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
641       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
642       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
643       ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr);
644       ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
645       /* Indent this deeper to emphasize the "inner" nature of this solver. */
646       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr);
647       ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr);
648       ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
649 
650       ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr);
651       if (sp) {
652         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
653       }
654 
655       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
656       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
657       if (flg) {
658         DM dmInner;
659 
660         /* Set DM for new solver */
661         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
662         ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr);
663         ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
664         ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr);
665         ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
666       } else {
667         ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr);
668         ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
669       }
670       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
671 
672       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
673       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
674       if (flg) {
675         DM dmInner;
676 
677         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
678         ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr);
679         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
680         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
681         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
682         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
683         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
684         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr);
685         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
686       } else {
687         jac->kspupper = jac->head->ksp;
688         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
689       }
690 
691       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
692       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
693       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
694       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
695       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
696         PC pcschur;
697         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
698         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
699         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
700       }
701       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
702       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
703       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
704       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
705       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
706     }
707 
708     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
709     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
710     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
711     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
712     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
713     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
714     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
715     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
716     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
717   } else {
718     /* set up the individual splits' PCs */
719     i     = 0;
720     ilink = jac->head;
721     while (ilink) {
722       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
723       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
724       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
725       i++;
726       ilink = ilink->next;
727     }
728   }
729 
730   jac->suboptionsset = PETSC_TRUE;
731   PetscFunctionReturn(0);
732 }
733 
734 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
735   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
736    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
737    KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
738    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
739    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
740 
741 #undef __FUNCT__
742 #define __FUNCT__ "PCApply_FieldSplit_Schur"
743 static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
744 {
745   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
746   PetscErrorCode    ierr;
747   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
748   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
749 
750   PetscFunctionBegin;
751   switch (jac->schurfactorization) {
752   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
753     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
754     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
755     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
756     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
757     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
758     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
759     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
760     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
761     ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
762     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
763     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
764     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
765     break;
766   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
767     /* [A00 0; A10 S], suitable for left preconditioning */
768     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
769     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
770     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
771     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
772     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
773     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
774     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
775     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
776     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
777     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
778     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
779     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
780     break;
781   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
782     /* [A00 A01; 0 S], suitable for right preconditioning */
783     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
784     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
785     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
786     ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
787     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
788     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
789     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
790     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
791     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
792     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
793     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
794     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
795     break;
796   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
797     /* [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 */
798     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
799     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
800     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
801     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
802     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
803     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
804     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
805 
806     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
807     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
808     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
809 
810     if (kspUpper == kspA) {
811       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
812       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
813       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
814     } else {
815       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
816       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
817       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
818       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
819     }
820     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
821     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
822   }
823   PetscFunctionReturn(0);
824 }
825 
826 #undef __FUNCT__
827 #define __FUNCT__ "PCApply_FieldSplit"
828 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
829 {
830   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
831   PetscErrorCode    ierr;
832   PC_FieldSplitLink ilink = jac->head;
833   PetscInt          cnt,bs;
834 
835   PetscFunctionBegin;
836   if (jac->type == PC_COMPOSITE_ADDITIVE) {
837     if (jac->defaultsplit) {
838       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
839       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);
840       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
841       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);
842       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
843       while (ilink) {
844         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
845         ilink = ilink->next;
846       }
847       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
848     } else {
849       ierr = VecSet(y,0.0);CHKERRQ(ierr);
850       while (ilink) {
851         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
852         ilink = ilink->next;
853       }
854     }
855   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
856     if (!jac->w1) {
857       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
858       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
859     }
860     ierr = VecSet(y,0.0);CHKERRQ(ierr);
861     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
862     cnt  = 1;
863     while (ilink->next) {
864       ilink = ilink->next;
865       /* compute the residual only over the part of the vector needed */
866       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
867       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
868       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
869       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
870       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
871       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
872       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
873     }
874     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
875       cnt -= 2;
876       while (ilink->previous) {
877         ilink = ilink->previous;
878         /* compute the residual only over the part of the vector needed */
879         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
880         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
881         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
882         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
883         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
884         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
885         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
886       }
887     }
888   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
889   PetscFunctionReturn(0);
890 }
891 
892 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
893   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
894    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
895    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
896    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
897    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
898 
899 #undef __FUNCT__
900 #define __FUNCT__ "PCApplyTranspose_FieldSplit"
901 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
902 {
903   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
904   PetscErrorCode    ierr;
905   PC_FieldSplitLink ilink = jac->head;
906   PetscInt          bs;
907 
908   PetscFunctionBegin;
909   if (jac->type == PC_COMPOSITE_ADDITIVE) {
910     if (jac->defaultsplit) {
911       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
912       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);
913       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
914       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);
915       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
916       while (ilink) {
917         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
918         ilink = ilink->next;
919       }
920       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
921     } else {
922       ierr = VecSet(y,0.0);CHKERRQ(ierr);
923       while (ilink) {
924         ierr  = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
925         ilink = ilink->next;
926       }
927     }
928   } else {
929     if (!jac->w1) {
930       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
931       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
932     }
933     ierr = VecSet(y,0.0);CHKERRQ(ierr);
934     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
935       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
936       while (ilink->next) {
937         ilink = ilink->next;
938         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
939         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
940         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
941       }
942       while (ilink->previous) {
943         ilink = ilink->previous;
944         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
945         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
946         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
947       }
948     } else {
949       while (ilink->next) {   /* get to last entry in linked list */
950         ilink = ilink->next;
951       }
952       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
953       while (ilink->previous) {
954         ilink = ilink->previous;
955         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
956         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
957         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
958       }
959     }
960   }
961   PetscFunctionReturn(0);
962 }
963 
964 #undef __FUNCT__
965 #define __FUNCT__ "PCReset_FieldSplit"
966 static PetscErrorCode PCReset_FieldSplit(PC pc)
967 {
968   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
969   PetscErrorCode    ierr;
970   PC_FieldSplitLink ilink = jac->head,next;
971 
972   PetscFunctionBegin;
973   while (ilink) {
974     ierr  = KSPReset(ilink->ksp);CHKERRQ(ierr);
975     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
976     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
977     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
978     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
979     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
980     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
981     next  = ilink->next;
982     ilink = next;
983   }
984   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
985   if (jac->mat && jac->mat != jac->pmat) {
986     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
987   } else if (jac->mat) {
988     jac->mat = NULL;
989   }
990   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
991   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
992   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
993   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
994   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
995   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
996   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
997   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
998   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
999   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
1000   jac->reset = PETSC_TRUE;
1001   PetscFunctionReturn(0);
1002 }
1003 
1004 #undef __FUNCT__
1005 #define __FUNCT__ "PCDestroy_FieldSplit"
1006 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1007 {
1008   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1009   PetscErrorCode    ierr;
1010   PC_FieldSplitLink ilink = jac->head,next;
1011 
1012   PetscFunctionBegin;
1013   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1014   while (ilink) {
1015     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1016     next  = ilink->next;
1017     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1018     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
1019     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
1020     ierr  = PetscFree(ilink);CHKERRQ(ierr);
1021     ilink = next;
1022   }
1023   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1024   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1025   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr);
1026   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr);
1027   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr);
1028   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr);
1030   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",NULL);CHKERRQ(ierr);
1031   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr);
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 #undef __FUNCT__
1036 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
1037 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
1038 {
1039   PetscErrorCode  ierr;
1040   PetscInt        bs;
1041   PetscBool       flg,stokes = PETSC_FALSE;
1042   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
1043   PCCompositeType ctype;
1044 
1045   PetscFunctionBegin;
1046   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
1047   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
1048   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
1049   if (flg) {
1050     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
1051   }
1052 
1053   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1054   if (stokes) {
1055     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1056     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1057   }
1058 
1059   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
1060   if (flg) {
1061     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
1062   }
1063   /* Only setup fields once */
1064   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1065     /* only allow user to set fields from command line if bs is already known.
1066        otherwise user can set them in PCFieldSplitSetDefaults() */
1067     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
1068     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1069   }
1070   if (jac->type == PC_COMPOSITE_SCHUR) {
1071     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1072     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
1073     ierr = 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);CHKERRQ(ierr);
1074     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);CHKERRQ(ierr);
1075   }
1076   ierr = PetscOptionsTail();CHKERRQ(ierr);
1077   PetscFunctionReturn(0);
1078 }
1079 
1080 /*------------------------------------------------------------------------------------*/
1081 
1082 #undef __FUNCT__
1083 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
1084 static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1085 {
1086   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1087   PetscErrorCode    ierr;
1088   PC_FieldSplitLink ilink,next = jac->head;
1089   char              prefix[128];
1090   PetscInt          i;
1091 
1092   PetscFunctionBegin;
1093   if (jac->splitdefined) {
1094     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
1095     PetscFunctionReturn(0);
1096   }
1097   for (i=0; i<n; i++) {
1098     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1099     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1100   }
1101   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1102   if (splitname) {
1103     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1104   } else {
1105     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1106     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1107   }
1108   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
1109   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1110   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
1111   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
1112 
1113   ilink->nfields = n;
1114   ilink->next    = NULL;
1115   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1116   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1117   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
1118   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1119 
1120   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1121   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1122 
1123   if (!next) {
1124     jac->head       = ilink;
1125     ilink->previous = NULL;
1126   } else {
1127     while (next->next) {
1128       next = next->next;
1129     }
1130     next->next      = ilink;
1131     ilink->previous = next;
1132   }
1133   jac->nsplits++;
1134   PetscFunctionReturn(0);
1135 }
1136 
1137 #undef __FUNCT__
1138 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
1139 static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1140 {
1141   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1142   PetscErrorCode ierr;
1143 
1144   PetscFunctionBegin;
1145   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1146   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1147 
1148   (*subksp)[1] = jac->kspschur;
1149   if (n) *n = jac->nsplits;
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
1155 static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1156 {
1157   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1158   PetscErrorCode    ierr;
1159   PetscInt          cnt   = 0;
1160   PC_FieldSplitLink ilink = jac->head;
1161 
1162   PetscFunctionBegin;
1163   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1164   while (ilink) {
1165     (*subksp)[cnt++] = ilink->ksp;
1166     ilink            = ilink->next;
1167   }
1168   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);
1169   if (n) *n = jac->nsplits;
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 #undef __FUNCT__
1174 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
1175 static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1176 {
1177   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1178   PetscErrorCode    ierr;
1179   PC_FieldSplitLink ilink, next = jac->head;
1180   char              prefix[128];
1181 
1182   PetscFunctionBegin;
1183   if (jac->splitdefined) {
1184     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
1185     PetscFunctionReturn(0);
1186   }
1187   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1188   if (splitname) {
1189     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1190   } else {
1191     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1192     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1193   }
1194   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1195   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1196   ilink->is     = is;
1197   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1198   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1199   ilink->is_col = is;
1200   ilink->next   = NULL;
1201   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1202   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1203   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
1204   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1205 
1206   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1207   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1208 
1209   if (!next) {
1210     jac->head       = ilink;
1211     ilink->previous = NULL;
1212   } else {
1213     while (next->next) {
1214       next = next->next;
1215     }
1216     next->next      = ilink;
1217     ilink->previous = next;
1218   }
1219   jac->nsplits++;
1220   PetscFunctionReturn(0);
1221 }
1222 
1223 #undef __FUNCT__
1224 #define __FUNCT__ "PCFieldSplitSetFields"
1225 /*@
1226     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
1227 
1228     Logically Collective on PC
1229 
1230     Input Parameters:
1231 +   pc  - the preconditioner context
1232 .   splitname - name of this split, if NULL the number of the split is used
1233 .   n - the number of fields in this split
1234 -   fields - the fields in this split
1235 
1236     Level: intermediate
1237 
1238     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1239 
1240      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1241      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
1242      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1243      where the numbered entries indicate what is in the field.
1244 
1245      This function is called once per split (it creates a new split each time).  Solve options
1246      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1247 
1248      Developer Note: This routine does not actually create the IS representing the split, that is delayed
1249      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
1250      available when this routine is called.
1251 
1252 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
1253 
1254 @*/
1255 PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1256 {
1257   PetscErrorCode ierr;
1258 
1259   PetscFunctionBegin;
1260   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1261   PetscValidCharPointer(splitname,2);
1262   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1263   PetscValidIntPointer(fields,3);
1264   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
1265   PetscFunctionReturn(0);
1266 }
1267 
1268 #undef __FUNCT__
1269 #define __FUNCT__ "PCFieldSplitSetIS"
1270 /*@
1271     PCFieldSplitSetIS - Sets the exact elements for field
1272 
1273     Logically Collective on PC
1274 
1275     Input Parameters:
1276 +   pc  - the preconditioner context
1277 .   splitname - name of this split, if NULL the number of the split is used
1278 -   is - the index set that defines the vector elements in this field
1279 
1280 
1281     Notes:
1282     Use PCFieldSplitSetFields(), for fields defined by strided types.
1283 
1284     This function is called once per split (it creates a new split each time).  Solve options
1285     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1286 
1287     Level: intermediate
1288 
1289 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1290 
1291 @*/
1292 PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1293 {
1294   PetscErrorCode ierr;
1295 
1296   PetscFunctionBegin;
1297   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1298   if (splitname) PetscValidCharPointer(splitname,2);
1299   PetscValidHeaderSpecific(is,IS_CLASSID,3);
1300   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 #undef __FUNCT__
1305 #define __FUNCT__ "PCFieldSplitGetIS"
1306 /*@
1307     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
1308 
1309     Logically Collective on PC
1310 
1311     Input Parameters:
1312 +   pc  - the preconditioner context
1313 -   splitname - name of this split
1314 
1315     Output Parameter:
1316 -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
1317 
1318     Level: intermediate
1319 
1320 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
1321 
1322 @*/
1323 PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
1324 {
1325   PetscErrorCode ierr;
1326 
1327   PetscFunctionBegin;
1328   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1329   PetscValidCharPointer(splitname,2);
1330   PetscValidPointer(is,3);
1331   {
1332     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
1333     PC_FieldSplitLink ilink = jac->head;
1334     PetscBool         found;
1335 
1336     *is = NULL;
1337     while (ilink) {
1338       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
1339       if (found) {
1340         *is = ilink->is;
1341         break;
1342       }
1343       ilink = ilink->next;
1344     }
1345   }
1346   PetscFunctionReturn(0);
1347 }
1348 
1349 #undef __FUNCT__
1350 #define __FUNCT__ "PCFieldSplitSetBlockSize"
1351 /*@
1352     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
1353       fieldsplit preconditioner. If not set the matrix block size is used.
1354 
1355     Logically Collective on PC
1356 
1357     Input Parameters:
1358 +   pc  - the preconditioner context
1359 -   bs - the block size
1360 
1361     Level: intermediate
1362 
1363 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
1364 
1365 @*/
1366 PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1367 {
1368   PetscErrorCode ierr;
1369 
1370   PetscFunctionBegin;
1371   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1372   PetscValidLogicalCollectiveInt(pc,bs,2);
1373   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
1374   PetscFunctionReturn(0);
1375 }
1376 
1377 #undef __FUNCT__
1378 #define __FUNCT__ "PCFieldSplitGetSubKSP"
1379 /*@C
1380    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
1381 
1382    Collective on KSP
1383 
1384    Input Parameter:
1385 .  pc - the preconditioner context
1386 
1387    Output Parameters:
1388 +  n - the number of splits
1389 -  pc - the array of KSP contexts
1390 
1391    Note:
1392    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1393    (not the KSP just the array that contains them).
1394 
1395    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
1396 
1397    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1398       You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the
1399       KSP array must be.
1400 
1401 
1402    Level: advanced
1403 
1404 .seealso: PCFIELDSPLIT
1405 @*/
1406 PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1407 {
1408   PetscErrorCode ierr;
1409 
1410   PetscFunctionBegin;
1411   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1412   if (n) PetscValidIntPointer(n,2);
1413   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
1414   PetscFunctionReturn(0);
1415 }
1416 
1417 #undef __FUNCT__
1418 #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1419 /*@
1420     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1421       A11 matrix. Otherwise no preconditioner is used.
1422 
1423     Collective on PC
1424 
1425     Input Parameters:
1426 +   pc  - the preconditioner context
1427 .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default
1428 -   userpre - matrix to use for preconditioning, or NULL
1429 
1430     Options Database:
1431 .     -pc_fieldsplit_schur_precondition <self,user,a11> default is a11
1432 
1433     Notes:
1434 $    If ptype is
1435 $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1436 $             to this function).
1437 $        a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1438 $             matrix associated with the Schur complement (i.e. A11)
1439 $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1440 $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1441 $             preconditioner
1442 
1443      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1444     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1445     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1446 
1447     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1448      the name since it sets a proceedure to use.
1449 
1450     Level: intermediate
1451 
1452 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1453 
1454 @*/
1455 PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1456 {
1457   PetscErrorCode ierr;
1458 
1459   PetscFunctionBegin;
1460   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1461   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1462   PetscFunctionReturn(0);
1463 }
1464 
1465 #undef __FUNCT__
1466 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1467 static PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1468 {
1469   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1470   PetscErrorCode ierr;
1471 
1472   PetscFunctionBegin;
1473   jac->schurpre = ptype;
1474   if (pre) {
1475     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1476     jac->schur_user = pre;
1477     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1478   }
1479   PetscFunctionReturn(0);
1480 }
1481 
1482 #undef __FUNCT__
1483 #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1484 /*@
1485     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1486 
1487     Collective on PC
1488 
1489     Input Parameters:
1490 +   pc  - the preconditioner context
1491 -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1492 
1493     Options Database:
1494 .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1495 
1496 
1497     Level: intermediate
1498 
1499     Notes:
1500     The FULL factorization is
1501 
1502 $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1503 $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1504 
1505     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,
1506     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).
1507 
1508     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
1509     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
1510     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,
1511     the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.
1512 
1513     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
1514     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).
1515 
1516     References:
1517     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1518 
1519     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1520 
1521 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1522 @*/
1523 PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1524 {
1525   PetscErrorCode ierr;
1526 
1527   PetscFunctionBegin;
1528   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1529   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1530   PetscFunctionReturn(0);
1531 }
1532 
1533 #undef __FUNCT__
1534 #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1535 static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1536 {
1537   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1538 
1539   PetscFunctionBegin;
1540   jac->schurfactorization = ftype;
1541   PetscFunctionReturn(0);
1542 }
1543 
1544 #undef __FUNCT__
1545 #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
1546 /*@C
1547    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
1548 
1549    Collective on KSP
1550 
1551    Input Parameter:
1552 .  pc - the preconditioner context
1553 
1554    Output Parameters:
1555 +  A00 - the (0,0) block
1556 .  A01 - the (0,1) block
1557 .  A10 - the (1,0) block
1558 -  A11 - the (1,1) block
1559 
1560    Level: advanced
1561 
1562 .seealso: PCFIELDSPLIT
1563 @*/
1564 PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
1565 {
1566   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
1567 
1568   PetscFunctionBegin;
1569   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1570   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1571   if (A00) *A00 = jac->pmat[0];
1572   if (A01) *A01 = jac->B;
1573   if (A10) *A10 = jac->C;
1574   if (A11) *A11 = jac->pmat[1];
1575   PetscFunctionReturn(0);
1576 }
1577 
1578 #undef __FUNCT__
1579 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1580 static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
1581 {
1582   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1583   PetscErrorCode ierr;
1584 
1585   PetscFunctionBegin;
1586   jac->type = type;
1587   if (type == PC_COMPOSITE_SCHUR) {
1588     pc->ops->apply = PCApply_FieldSplit_Schur;
1589     pc->ops->view  = PCView_FieldSplit_Schur;
1590 
1591     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1592     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1593     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1594 
1595   } else {
1596     pc->ops->apply = PCApply_FieldSplit;
1597     pc->ops->view  = PCView_FieldSplit;
1598 
1599     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1600     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",0);CHKERRQ(ierr);
1601     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
1602   }
1603   PetscFunctionReturn(0);
1604 }
1605 
1606 #undef __FUNCT__
1607 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
1608 static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
1609 {
1610   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1611 
1612   PetscFunctionBegin;
1613   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1614   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);
1615   jac->bs = bs;
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 #undef __FUNCT__
1620 #define __FUNCT__ "PCFieldSplitSetType"
1621 /*@
1622    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
1623 
1624    Collective on PC
1625 
1626    Input Parameter:
1627 .  pc - the preconditioner context
1628 .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1629 
1630    Options Database Key:
1631 .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
1632 
1633    Level: Intermediate
1634 
1635 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1636 
1637 .seealso: PCCompositeSetType()
1638 
1639 @*/
1640 PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
1641 {
1642   PetscErrorCode ierr;
1643 
1644   PetscFunctionBegin;
1645   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1646   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
1647   PetscFunctionReturn(0);
1648 }
1649 
1650 #undef __FUNCT__
1651 #define __FUNCT__ "PCFieldSplitGetType"
1652 /*@
1653   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
1654 
1655   Not collective
1656 
1657   Input Parameter:
1658 . pc - the preconditioner context
1659 
1660   Output Parameter:
1661 . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1662 
1663   Level: Intermediate
1664 
1665 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1666 .seealso: PCCompositeSetType()
1667 @*/
1668 PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
1669 {
1670   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
1671 
1672   PetscFunctionBegin;
1673   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1674   PetscValidIntPointer(type,2);
1675   *type = jac->type;
1676   PetscFunctionReturn(0);
1677 }
1678 
1679 #undef __FUNCT__
1680 #define __FUNCT__ "PCFieldSplitSetDMSplits"
1681 /*@
1682    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
1683 
1684    Logically Collective
1685 
1686    Input Parameters:
1687 +  pc   - the preconditioner context
1688 -  flg  - boolean indicating whether to use field splits defined by the DM
1689 
1690    Options Database Key:
1691 .  -pc_fieldsplit_dm_splits
1692 
1693    Level: Intermediate
1694 
1695 .keywords: PC, DM, composite preconditioner, additive, multiplicative
1696 
1697 .seealso: PCFieldSplitGetDMSplits()
1698 
1699 @*/
1700 PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
1701 {
1702   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1703   PetscBool      isfs;
1704   PetscErrorCode ierr;
1705 
1706   PetscFunctionBegin;
1707   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1708   PetscValidLogicalCollectiveBool(pc,flg,2);
1709   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1710   if (isfs) {
1711     jac->dm_splits = flg;
1712   }
1713   PetscFunctionReturn(0);
1714 }
1715 
1716 
1717 #undef __FUNCT__
1718 #define __FUNCT__ "PCFieldSplitGetDMSplits"
1719 /*@
1720    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
1721 
1722    Logically Collective
1723 
1724    Input Parameter:
1725 .  pc   - the preconditioner context
1726 
1727    Output Parameter:
1728 .  flg  - boolean indicating whether to use field splits defined by the DM
1729 
1730    Level: Intermediate
1731 
1732 .keywords: PC, DM, composite preconditioner, additive, multiplicative
1733 
1734 .seealso: PCFieldSplitSetDMSplits()
1735 
1736 @*/
1737 PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
1738 {
1739   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1740   PetscBool      isfs;
1741   PetscErrorCode ierr;
1742 
1743   PetscFunctionBegin;
1744   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1745   PetscValidPointer(flg,2);
1746   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1747   if (isfs) {
1748     if(flg) *flg = jac->dm_splits;
1749   }
1750   PetscFunctionReturn(0);
1751 }
1752 
1753 /* -------------------------------------------------------------------------------------*/
1754 /*MC
1755    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1756                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
1757 
1758      To set options on the solvers for each block append -fieldsplit_ to all the PC
1759         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
1760 
1761      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
1762          and set the options directly on the resulting KSP object
1763 
1764    Level: intermediate
1765 
1766    Options Database Keys:
1767 +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
1768 .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
1769                               been supplied explicitly by -pc_fieldsplit_%d_fields
1770 .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
1771 .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
1772 .   -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11
1773 .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
1774 
1775 -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
1776      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
1777 
1778    Notes:
1779     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1780      to define a field by an arbitrary collection of entries.
1781 
1782       If no fields are set the default is used. The fields are defined by entries strided by bs,
1783       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1784       if this is not called the block size defaults to the blocksize of the second matrix passed
1785       to KSPSetOperators()/PCSetOperators().
1786 
1787 $     For the Schur complement preconditioner if J = ( A00 A01 )
1788 $                                                    ( A10 A11 )
1789 $     the preconditioner using full factorization is
1790 $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1791 $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1792      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1793      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1794      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1795      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1796      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1797      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1798      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
1799      diag gives
1800 $              ( inv(A00)     0   )
1801 $              (   0      -ksp(S) )
1802      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
1803 $              (  A00   0 )
1804 $              (  A10   S )
1805      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1806 $              ( A00 A01 )
1807 $              (  0   S  )
1808      where again the inverses of A00 and S are applied using KSPs.
1809 
1810      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1811      is used automatically for a second block.
1812 
1813      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1814      Generally it should be used with the AIJ format.
1815 
1816      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1817      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1818      inside a smoother resulting in "Distributive Smoothers".
1819 
1820    Concepts: physics based preconditioners, block preconditioners
1821 
1822 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1823            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1824 M*/
1825 
1826 #undef __FUNCT__
1827 #define __FUNCT__ "PCCreate_FieldSplit"
1828 PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
1829 {
1830   PetscErrorCode ierr;
1831   PC_FieldSplit  *jac;
1832 
1833   PetscFunctionBegin;
1834   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
1835 
1836   jac->bs                 = -1;
1837   jac->nsplits            = 0;
1838   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
1839   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1840   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
1841   jac->dm_splits          = PETSC_TRUE;
1842 
1843   pc->data = (void*)jac;
1844 
1845   pc->ops->apply           = PCApply_FieldSplit;
1846   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
1847   pc->ops->setup           = PCSetUp_FieldSplit;
1848   pc->ops->reset           = PCReset_FieldSplit;
1849   pc->ops->destroy         = PCDestroy_FieldSplit;
1850   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
1851   pc->ops->view            = PCView_FieldSplit;
1852   pc->ops->applyrichardson = 0;
1853 
1854   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1855   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1856   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
1857   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
1858   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
1859   PetscFunctionReturn(0);
1860 }
1861 
1862 
1863 
1864