xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision eb975fb2d04d974246739c9c552b1cfae0c3a9e4)
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   CHKMEMQ;
837   if (jac->type == PC_COMPOSITE_ADDITIVE) {
838     if (jac->defaultsplit) {
839       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
840       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);
841       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
842       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);
843       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
844       while (ilink) {
845         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
846         ilink = ilink->next;
847       }
848       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
849     } else {
850       ierr = VecSet(y,0.0);CHKERRQ(ierr);
851       while (ilink) {
852         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
853         ilink = ilink->next;
854       }
855     }
856   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
857     if (!jac->w1) {
858       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
859       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
860     }
861     ierr = VecSet(y,0.0);CHKERRQ(ierr);
862     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
863     cnt  = 1;
864     while (ilink->next) {
865       ilink = ilink->next;
866       /* compute the residual only over the part of the vector needed */
867       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
868       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
869       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
870       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
871       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
872       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
873       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
874     }
875     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
876       cnt -= 2;
877       while (ilink->previous) {
878         ilink = ilink->previous;
879         /* compute the residual only over the part of the vector needed */
880         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
881         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
882         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
883         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
884         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
885         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
886         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
887       }
888     }
889   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
890   CHKMEMQ;
891   PetscFunctionReturn(0);
892 }
893 
894 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
895   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
896    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
897    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
898    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
899    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
900 
901 #undef __FUNCT__
902 #define __FUNCT__ "PCApplyTranspose_FieldSplit"
903 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
904 {
905   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
906   PetscErrorCode    ierr;
907   PC_FieldSplitLink ilink = jac->head;
908   PetscInt          bs;
909 
910   PetscFunctionBegin;
911   CHKMEMQ;
912   if (jac->type == PC_COMPOSITE_ADDITIVE) {
913     if (jac->defaultsplit) {
914       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
915       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);
916       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
917       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);
918       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
919       while (ilink) {
920         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
921         ilink = ilink->next;
922       }
923       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
924     } else {
925       ierr = VecSet(y,0.0);CHKERRQ(ierr);
926       while (ilink) {
927         ierr  = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
928         ilink = ilink->next;
929       }
930     }
931   } else {
932     if (!jac->w1) {
933       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
934       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
935     }
936     ierr = VecSet(y,0.0);CHKERRQ(ierr);
937     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
938       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
939       while (ilink->next) {
940         ilink = ilink->next;
941         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
942         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
943         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
944       }
945       while (ilink->previous) {
946         ilink = ilink->previous;
947         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
948         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
949         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
950       }
951     } else {
952       while (ilink->next) {   /* get to last entry in linked list */
953         ilink = ilink->next;
954       }
955       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
956       while (ilink->previous) {
957         ilink = ilink->previous;
958         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
959         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
960         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
961       }
962     }
963   }
964   CHKMEMQ;
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNCT__
969 #define __FUNCT__ "PCReset_FieldSplit"
970 static PetscErrorCode PCReset_FieldSplit(PC pc)
971 {
972   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
973   PetscErrorCode    ierr;
974   PC_FieldSplitLink ilink = jac->head,next;
975 
976   PetscFunctionBegin;
977   while (ilink) {
978     ierr  = KSPReset(ilink->ksp);CHKERRQ(ierr);
979     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
980     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
981     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
982     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
983     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
984     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
985     next  = ilink->next;
986     ilink = next;
987   }
988   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
989   if (jac->mat && jac->mat != jac->pmat) {
990     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
991   } else if (jac->mat) {
992     jac->mat = NULL;
993   }
994   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
995   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
996   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
997   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
998   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
999   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1000   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1001   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
1002   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
1003   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
1004   jac->reset = PETSC_TRUE;
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNCT__
1009 #define __FUNCT__ "PCDestroy_FieldSplit"
1010 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1011 {
1012   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1013   PetscErrorCode    ierr;
1014   PC_FieldSplitLink ilink = jac->head,next;
1015 
1016   PetscFunctionBegin;
1017   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1018   while (ilink) {
1019     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1020     next  = ilink->next;
1021     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1022     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
1023     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
1024     ierr  = PetscFree(ilink);CHKERRQ(ierr);
1025     ilink = next;
1026   }
1027   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1028   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",NULL);CHKERRQ(ierr);
1030   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C","",NULL);CHKERRQ(ierr);
1031   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C","",NULL);CHKERRQ(ierr);
1032   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C","",NULL);CHKERRQ(ierr);
1033   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",NULL);CHKERRQ(ierr);
1034   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",NULL);CHKERRQ(ierr);
1035   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",NULL);CHKERRQ(ierr);
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 #undef __FUNCT__
1040 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
1041 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
1042 {
1043   PetscErrorCode  ierr;
1044   PetscInt        bs;
1045   PetscBool       flg,stokes = PETSC_FALSE;
1046   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
1047   PCCompositeType ctype;
1048 
1049   PetscFunctionBegin;
1050   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
1051   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
1052   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
1053   if (flg) {
1054     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
1055   }
1056 
1057   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1058   if (stokes) {
1059     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1060     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1061   }
1062 
1063   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
1064   if (flg) {
1065     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
1066   }
1067   /* Only setup fields once */
1068   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1069     /* only allow user to set fields from command line if bs is already known.
1070        otherwise user can set them in PCFieldSplitSetDefaults() */
1071     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
1072     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1073   }
1074   if (jac->type == PC_COMPOSITE_SCHUR) {
1075     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1076     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
1077     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);
1078     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);CHKERRQ(ierr);
1079   }
1080   ierr = PetscOptionsTail();CHKERRQ(ierr);
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 /*------------------------------------------------------------------------------------*/
1085 
1086 #undef __FUNCT__
1087 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
1088 static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1089 {
1090   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1091   PetscErrorCode    ierr;
1092   PC_FieldSplitLink ilink,next = jac->head;
1093   char              prefix[128];
1094   PetscInt          i;
1095 
1096   PetscFunctionBegin;
1097   if (jac->splitdefined) {
1098     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
1099     PetscFunctionReturn(0);
1100   }
1101   for (i=0; i<n; i++) {
1102     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1103     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1104   }
1105   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1106   if (splitname) {
1107     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1108   } else {
1109     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1110     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1111   }
1112   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
1113   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1114   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
1115   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
1116 
1117   ilink->nfields = n;
1118   ilink->next    = NULL;
1119   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1120   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1121   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
1122   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1123 
1124   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1125   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1126 
1127   if (!next) {
1128     jac->head       = ilink;
1129     ilink->previous = NULL;
1130   } else {
1131     while (next->next) {
1132       next = next->next;
1133     }
1134     next->next      = ilink;
1135     ilink->previous = next;
1136   }
1137   jac->nsplits++;
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 #undef __FUNCT__
1142 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
1143 static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1144 {
1145   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1146   PetscErrorCode ierr;
1147 
1148   PetscFunctionBegin;
1149   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1150   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1151 
1152   (*subksp)[1] = jac->kspschur;
1153   if (n) *n = jac->nsplits;
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 #undef __FUNCT__
1158 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
1159 static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1160 {
1161   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1162   PetscErrorCode    ierr;
1163   PetscInt          cnt   = 0;
1164   PC_FieldSplitLink ilink = jac->head;
1165 
1166   PetscFunctionBegin;
1167   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1168   while (ilink) {
1169     (*subksp)[cnt++] = ilink->ksp;
1170     ilink            = ilink->next;
1171   }
1172   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);
1173   if (n) *n = jac->nsplits;
1174   PetscFunctionReturn(0);
1175 }
1176 
1177 #undef __FUNCT__
1178 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
1179 static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1180 {
1181   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1182   PetscErrorCode    ierr;
1183   PC_FieldSplitLink ilink, next = jac->head;
1184   char              prefix[128];
1185 
1186   PetscFunctionBegin;
1187   if (jac->splitdefined) {
1188     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
1189     PetscFunctionReturn(0);
1190   }
1191   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1192   if (splitname) {
1193     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1194   } else {
1195     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1196     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1197   }
1198   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1199   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1200   ilink->is     = is;
1201   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1202   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1203   ilink->is_col = is;
1204   ilink->next   = NULL;
1205   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1206   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1207   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
1208   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1209 
1210   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1211   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1212 
1213   if (!next) {
1214     jac->head       = ilink;
1215     ilink->previous = NULL;
1216   } else {
1217     while (next->next) {
1218       next = next->next;
1219     }
1220     next->next      = ilink;
1221     ilink->previous = next;
1222   }
1223   jac->nsplits++;
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 #undef __FUNCT__
1228 #define __FUNCT__ "PCFieldSplitSetFields"
1229 /*@
1230     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
1231 
1232     Logically Collective on PC
1233 
1234     Input Parameters:
1235 +   pc  - the preconditioner context
1236 .   splitname - name of this split, if NULL the number of the split is used
1237 .   n - the number of fields in this split
1238 -   fields - the fields in this split
1239 
1240     Level: intermediate
1241 
1242     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1243 
1244      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1245      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
1246      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1247      where the numbered entries indicate what is in the field.
1248 
1249      This function is called once per split (it creates a new split each time).  Solve options
1250      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1251 
1252      Developer Note: This routine does not actually create the IS representing the split, that is delayed
1253      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
1254      available when this routine is called.
1255 
1256 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
1257 
1258 @*/
1259 PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1260 {
1261   PetscErrorCode ierr;
1262 
1263   PetscFunctionBegin;
1264   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1265   PetscValidCharPointer(splitname,2);
1266   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1267   PetscValidIntPointer(fields,3);
1268   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 #undef __FUNCT__
1273 #define __FUNCT__ "PCFieldSplitSetIS"
1274 /*@
1275     PCFieldSplitSetIS - Sets the exact elements for field
1276 
1277     Logically Collective on PC
1278 
1279     Input Parameters:
1280 +   pc  - the preconditioner context
1281 .   splitname - name of this split, if NULL the number of the split is used
1282 -   is - the index set that defines the vector elements in this field
1283 
1284 
1285     Notes:
1286     Use PCFieldSplitSetFields(), for fields defined by strided types.
1287 
1288     This function is called once per split (it creates a new split each time).  Solve options
1289     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1290 
1291     Level: intermediate
1292 
1293 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1294 
1295 @*/
1296 PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1297 {
1298   PetscErrorCode ierr;
1299 
1300   PetscFunctionBegin;
1301   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1302   if (splitname) PetscValidCharPointer(splitname,2);
1303   PetscValidHeaderSpecific(is,IS_CLASSID,3);
1304   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 #undef __FUNCT__
1309 #define __FUNCT__ "PCFieldSplitGetIS"
1310 /*@
1311     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
1312 
1313     Logically Collective on PC
1314 
1315     Input Parameters:
1316 +   pc  - the preconditioner context
1317 -   splitname - name of this split
1318 
1319     Output Parameter:
1320 -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
1321 
1322     Level: intermediate
1323 
1324 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
1325 
1326 @*/
1327 PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
1328 {
1329   PetscErrorCode ierr;
1330 
1331   PetscFunctionBegin;
1332   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1333   PetscValidCharPointer(splitname,2);
1334   PetscValidPointer(is,3);
1335   {
1336     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
1337     PC_FieldSplitLink ilink = jac->head;
1338     PetscBool         found;
1339 
1340     *is = NULL;
1341     while (ilink) {
1342       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
1343       if (found) {
1344         *is = ilink->is;
1345         break;
1346       }
1347       ilink = ilink->next;
1348     }
1349   }
1350   PetscFunctionReturn(0);
1351 }
1352 
1353 #undef __FUNCT__
1354 #define __FUNCT__ "PCFieldSplitSetBlockSize"
1355 /*@
1356     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
1357       fieldsplit preconditioner. If not set the matrix block size is used.
1358 
1359     Logically Collective on PC
1360 
1361     Input Parameters:
1362 +   pc  - the preconditioner context
1363 -   bs - the block size
1364 
1365     Level: intermediate
1366 
1367 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
1368 
1369 @*/
1370 PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
1371 {
1372   PetscErrorCode ierr;
1373 
1374   PetscFunctionBegin;
1375   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1376   PetscValidLogicalCollectiveInt(pc,bs,2);
1377   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
1378   PetscFunctionReturn(0);
1379 }
1380 
1381 #undef __FUNCT__
1382 #define __FUNCT__ "PCFieldSplitGetSubKSP"
1383 /*@C
1384    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
1385 
1386    Collective on KSP
1387 
1388    Input Parameter:
1389 .  pc - the preconditioner context
1390 
1391    Output Parameters:
1392 +  n - the number of splits
1393 -  pc - the array of KSP contexts
1394 
1395    Note:
1396    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1397    (not the KSP just the array that contains them).
1398 
1399    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
1400 
1401    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1402       You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the
1403       KSP array must be.
1404 
1405 
1406    Level: advanced
1407 
1408 .seealso: PCFIELDSPLIT
1409 @*/
1410 PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
1411 {
1412   PetscErrorCode ierr;
1413 
1414   PetscFunctionBegin;
1415   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1416   if (n) PetscValidIntPointer(n,2);
1417   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
1418   PetscFunctionReturn(0);
1419 }
1420 
1421 #undef __FUNCT__
1422 #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1423 /*@
1424     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1425       A11 matrix. Otherwise no preconditioner is used.
1426 
1427     Collective on PC
1428 
1429     Input Parameters:
1430 +   pc  - the preconditioner context
1431 .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default
1432 -   userpre - matrix to use for preconditioning, or NULL
1433 
1434     Options Database:
1435 .     -pc_fieldsplit_schur_precondition <self,user,a11> default is a11
1436 
1437     Notes:
1438 $    If ptype is
1439 $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1440 $             to this function).
1441 $        a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1442 $             matrix associated with the Schur complement (i.e. A11)
1443 $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1444 $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1445 $             preconditioner
1446 
1447      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1448     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1449     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1450 
1451     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1452      the name since it sets a proceedure to use.
1453 
1454     Level: intermediate
1455 
1456 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1457 
1458 @*/
1459 PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1460 {
1461   PetscErrorCode ierr;
1462 
1463   PetscFunctionBegin;
1464   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1465   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1466   PetscFunctionReturn(0);
1467 }
1468 
1469 #undef __FUNCT__
1470 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1471 static PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1472 {
1473   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1474   PetscErrorCode ierr;
1475 
1476   PetscFunctionBegin;
1477   jac->schurpre = ptype;
1478   if (pre) {
1479     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1480     jac->schur_user = pre;
1481     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1482   }
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 #undef __FUNCT__
1487 #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1488 /*@
1489     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1490 
1491     Collective on PC
1492 
1493     Input Parameters:
1494 +   pc  - the preconditioner context
1495 -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1496 
1497     Options Database:
1498 .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1499 
1500 
1501     Level: intermediate
1502 
1503     Notes:
1504     The FULL factorization is
1505 
1506 $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1507 $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1508 
1509     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,
1510     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).
1511 
1512     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
1513     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
1514     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,
1515     the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.
1516 
1517     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
1518     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).
1519 
1520     References:
1521     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1522 
1523     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1524 
1525 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1526 @*/
1527 PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1528 {
1529   PetscErrorCode ierr;
1530 
1531   PetscFunctionBegin;
1532   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1533   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1534   PetscFunctionReturn(0);
1535 }
1536 
1537 #undef __FUNCT__
1538 #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1539 static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1540 {
1541   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1542 
1543   PetscFunctionBegin;
1544   jac->schurfactorization = ftype;
1545   PetscFunctionReturn(0);
1546 }
1547 
1548 #undef __FUNCT__
1549 #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
1550 /*@C
1551    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
1552 
1553    Collective on KSP
1554 
1555    Input Parameter:
1556 .  pc - the preconditioner context
1557 
1558    Output Parameters:
1559 +  A00 - the (0,0) block
1560 .  A01 - the (0,1) block
1561 .  A10 - the (1,0) block
1562 -  A11 - the (1,1) block
1563 
1564    Level: advanced
1565 
1566 .seealso: PCFIELDSPLIT
1567 @*/
1568 PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
1569 {
1570   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
1571 
1572   PetscFunctionBegin;
1573   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1574   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1575   if (A00) *A00 = jac->pmat[0];
1576   if (A01) *A01 = jac->B;
1577   if (A10) *A10 = jac->C;
1578   if (A11) *A11 = jac->pmat[1];
1579   PetscFunctionReturn(0);
1580 }
1581 
1582 #undef __FUNCT__
1583 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1584 static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
1585 {
1586   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1587   PetscErrorCode ierr;
1588 
1589   PetscFunctionBegin;
1590   jac->type = type;
1591   if (type == PC_COMPOSITE_SCHUR) {
1592     pc->ops->apply = PCApply_FieldSplit_Schur;
1593     pc->ops->view  = PCView_FieldSplit_Schur;
1594 
1595     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1596     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1597     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1598 
1599   } else {
1600     pc->ops->apply = PCApply_FieldSplit;
1601     pc->ops->view  = PCView_FieldSplit;
1602 
1603     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1604     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1605     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
1606   }
1607   PetscFunctionReturn(0);
1608 }
1609 
1610 #undef __FUNCT__
1611 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
1612 static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
1613 {
1614   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1615 
1616   PetscFunctionBegin;
1617   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1618   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);
1619   jac->bs = bs;
1620   PetscFunctionReturn(0);
1621 }
1622 
1623 #undef __FUNCT__
1624 #define __FUNCT__ "PCFieldSplitSetType"
1625 /*@
1626    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
1627 
1628    Collective on PC
1629 
1630    Input Parameter:
1631 .  pc - the preconditioner context
1632 .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1633 
1634    Options Database Key:
1635 .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
1636 
1637    Level: Intermediate
1638 
1639 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1640 
1641 .seealso: PCCompositeSetType()
1642 
1643 @*/
1644 PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
1645 {
1646   PetscErrorCode ierr;
1647 
1648   PetscFunctionBegin;
1649   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1650   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
1651   PetscFunctionReturn(0);
1652 }
1653 
1654 #undef __FUNCT__
1655 #define __FUNCT__ "PCFieldSplitGetType"
1656 /*@
1657   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
1658 
1659   Not collective
1660 
1661   Input Parameter:
1662 . pc - the preconditioner context
1663 
1664   Output Parameter:
1665 . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1666 
1667   Level: Intermediate
1668 
1669 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1670 .seealso: PCCompositeSetType()
1671 @*/
1672 PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
1673 {
1674   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
1675 
1676   PetscFunctionBegin;
1677   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1678   PetscValidIntPointer(type,2);
1679   *type = jac->type;
1680   PetscFunctionReturn(0);
1681 }
1682 
1683 #undef __FUNCT__
1684 #define __FUNCT__ "PCFieldSplitSetDMSplits"
1685 /*@
1686    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
1687 
1688    Logically Collective
1689 
1690    Input Parameters:
1691 +  pc   - the preconditioner context
1692 -  flg  - boolean indicating whether to use field splits defined by the DM
1693 
1694    Options Database Key:
1695 .  -pc_fieldsplit_dm_splits
1696 
1697    Level: Intermediate
1698 
1699 .keywords: PC, DM, composite preconditioner, additive, multiplicative
1700 
1701 .seealso: PCFieldSplitGetDMSplits()
1702 
1703 @*/
1704 PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
1705 {
1706   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1707   PetscBool      isfs;
1708   PetscErrorCode ierr;
1709 
1710   PetscFunctionBegin;
1711   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1712   PetscValidLogicalCollectiveBool(pc,flg,2);
1713   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1714   if (isfs) {
1715     jac->dm_splits = flg;
1716   }
1717   PetscFunctionReturn(0);
1718 }
1719 
1720 
1721 #undef __FUNCT__
1722 #define __FUNCT__ "PCFieldSplitGetDMSplits"
1723 /*@
1724    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
1725 
1726    Logically Collective
1727 
1728    Input Parameter:
1729 .  pc   - the preconditioner context
1730 
1731    Output Parameter:
1732 .  flg  - boolean indicating whether to use field splits defined by the DM
1733 
1734    Level: Intermediate
1735 
1736 .keywords: PC, DM, composite preconditioner, additive, multiplicative
1737 
1738 .seealso: PCFieldSplitSetDMSplits()
1739 
1740 @*/
1741 PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
1742 {
1743   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1744   PetscBool      isfs;
1745   PetscErrorCode ierr;
1746 
1747   PetscFunctionBegin;
1748   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1749   PetscValidPointer(flg,2);
1750   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1751   if (isfs) {
1752     if(flg) *flg = jac->dm_splits;
1753   }
1754   PetscFunctionReturn(0);
1755 }
1756 
1757 /* -------------------------------------------------------------------------------------*/
1758 /*MC
1759    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1760                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
1761 
1762      To set options on the solvers for each block append -fieldsplit_ to all the PC
1763         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
1764 
1765      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
1766          and set the options directly on the resulting KSP object
1767 
1768    Level: intermediate
1769 
1770    Options Database Keys:
1771 +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
1772 .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
1773                               been supplied explicitly by -pc_fieldsplit_%d_fields
1774 .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
1775 .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
1776 .   -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11
1777 .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
1778 
1779 -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
1780      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
1781 
1782    Notes:
1783     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1784      to define a field by an arbitrary collection of entries.
1785 
1786       If no fields are set the default is used. The fields are defined by entries strided by bs,
1787       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1788       if this is not called the block size defaults to the blocksize of the second matrix passed
1789       to KSPSetOperators()/PCSetOperators().
1790 
1791 $     For the Schur complement preconditioner if J = ( A00 A01 )
1792 $                                                    ( A10 A11 )
1793 $     the preconditioner using full factorization is
1794 $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1795 $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1796      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1797      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1798      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1799      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1800      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1801      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1802      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
1803      diag gives
1804 $              ( inv(A00)     0   )
1805 $              (   0      -ksp(S) )
1806      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
1807 $              (  A00   0 )
1808 $              (  A10   S )
1809      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1810 $              ( A00 A01 )
1811 $              (  0   S  )
1812      where again the inverses of A00 and S are applied using KSPs.
1813 
1814      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1815      is used automatically for a second block.
1816 
1817      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1818      Generally it should be used with the AIJ format.
1819 
1820      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1821      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1822      inside a smoother resulting in "Distributive Smoothers".
1823 
1824    Concepts: physics based preconditioners, block preconditioners
1825 
1826 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1827            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1828 M*/
1829 
1830 #undef __FUNCT__
1831 #define __FUNCT__ "PCCreate_FieldSplit"
1832 PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
1833 {
1834   PetscErrorCode ierr;
1835   PC_FieldSplit  *jac;
1836 
1837   PetscFunctionBegin;
1838   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
1839 
1840   jac->bs                 = -1;
1841   jac->nsplits            = 0;
1842   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
1843   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1844   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
1845   jac->dm_splits          = PETSC_TRUE;
1846 
1847   pc->data = (void*)jac;
1848 
1849   pc->ops->apply           = PCApply_FieldSplit;
1850   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
1851   pc->ops->setup           = PCSetUp_FieldSplit;
1852   pc->ops->reset           = PCReset_FieldSplit;
1853   pc->ops->destroy         = PCDestroy_FieldSplit;
1854   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
1855   pc->ops->view            = PCView_FieldSplit;
1856   pc->ops->applyrichardson = 0;
1857 
1858   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1859   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1860   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
1861   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
1862   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
1863   PetscFunctionReturn(0);
1864 }
1865 
1866 
1867 
1868