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