xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 3e73264c0fb4e6a394f4588f9de8e8d942f17ae5)
1 #define PETSCKSP_DLL
2 
3 /*
4 
5 */
6 #include "private/pcimpl.h"     /*I "petscpc.h" I*/
7 
8 const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
9 
10 typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
11 struct _PC_FieldSplitLink {
12   KSP               ksp;
13   Vec               x,y;
14   PetscInt          nfields;
15   PetscInt          *fields;
16   VecScatter        sctx;
17   IS                is;
18   PC_FieldSplitLink next,previous;
19 };
20 
21 typedef struct {
22   PCCompositeType   type;
23   PetscTruth        defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
24   PetscTruth        realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
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   PetscTruth        issetup;
32   /* Only used when Schur complement preconditioning is used */
33   Mat               B;            /* The (0,1) block */
34   Mat               C;            /* The (1,0) block */
35   Mat               schur;        /* The Schur complement S = D - C A^{-1} B */
36   Mat               schur_user;   /* User-provided preconditioning matrix for the Schur complement */
37   PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
38   KSP               kspschur;     /* The solver for S */
39   PC_FieldSplitLink head;
40 } PC_FieldSplit;
41 
42 /*
43     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
44    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
45    PC you could change this.
46 */
47 
48 /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
49 * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
50 static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
51 {
52   switch (jac->schurpre) {
53     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
54     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
55     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
56     default:
57       return jac->schur_user ? jac->schur_user : jac->pmat[1];
58   }
59 }
60 
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "PCView_FieldSplit"
64 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
65 {
66   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
67   PetscErrorCode    ierr;
68   PetscTruth        iascii;
69   PetscInt          i,j;
70   PC_FieldSplitLink ilink = jac->head;
71 
72   PetscFunctionBegin;
73   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
74   if (iascii) {
75     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
76     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
77     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
78     for (i=0; i<jac->nsplits; i++) {
79       if (ilink->fields) {
80 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
81 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
82 	for (j=0; j<ilink->nfields; j++) {
83 	  if (j > 0) {
84 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
85 	  }
86 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
87 	}
88 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
89         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
90       } else {
91 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
92       }
93       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
94       ilink = ilink->next;
95     }
96     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
97   } else {
98     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 #undef __FUNCT__
104 #define __FUNCT__ "PCView_FieldSplit_Schur"
105 static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
106 {
107   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
108   PetscErrorCode    ierr;
109   PetscTruth        iascii;
110   PetscInt          i,j;
111   PC_FieldSplitLink ilink = jac->head;
112   KSP               ksp;
113 
114   PetscFunctionBegin;
115   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
116   if (iascii) {
117     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr);
118     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
119     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
120     for (i=0; i<jac->nsplits; i++) {
121       if (ilink->fields) {
122 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
123 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
124 	for (j=0; j<ilink->nfields; j++) {
125 	  if (j > 0) {
126 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
127 	  }
128 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
129 	}
130 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
131         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
132       } else {
133 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
134       }
135       ilink = ilink->next;
136     }
137     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
138     ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
139     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
140     ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
141     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
142     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
143     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
144     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
145     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
146     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
147   } else {
148     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
149   }
150   PetscFunctionReturn(0);
151 }
152 
153 #undef __FUNCT__
154 #define __FUNCT__ "PCFieldSplitSetDefaults"
155 static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
156 {
157   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
158   PetscErrorCode    ierr;
159   PC_FieldSplitLink ilink = jac->head;
160   PetscInt          i = 0,*ifields,nfields;
161   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
162   char              optionname[128];
163 
164   PetscFunctionBegin;
165   if (!ilink) {
166 
167     if (jac->bs <= 0) {
168       if (pc->pmat) {
169         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
170       } else {
171         jac->bs = 1;
172       }
173     }
174 
175     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
176     if (!flg) {
177       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
178          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
179       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
180       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
181       while (PETSC_TRUE) {
182         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
183         nfields = jac->bs;
184         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
185         if (!flg2) break;
186         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
187         flg = PETSC_FALSE;
188         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
189       }
190       ierr = PetscFree(ifields);CHKERRQ(ierr);
191     }
192 
193     if (flg) {
194       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
195       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
196       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
197       while (ilink) {
198 	for (i=0; i<ilink->nfields; i++) {
199 	  fields[ilink->fields[i]] = PETSC_TRUE;
200 	}
201 	ilink = ilink->next;
202       }
203       jac->defaultsplit = PETSC_TRUE;
204       for (i=0; i<jac->bs; i++) {
205 	if (!fields[i]) {
206 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
207 	} else {
208 	  jac->defaultsplit = PETSC_FALSE;
209 	}
210       }
211       ierr = PetscFree(fields);CHKERRQ(ierr);
212     }
213   } else if (jac->nsplits == 1) {
214     if (ilink->is) {
215       IS       is2;
216       PetscInt nmin,nmax;
217 
218       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
219       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
220       ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr);
221       ierr = ISDestroy(is2);CHKERRQ(ierr);
222     } else {
223       SETERRQ(PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
224     }
225   }
226   PetscFunctionReturn(0);
227 }
228 
229 
230 #undef __FUNCT__
231 #define __FUNCT__ "PCSetUp_FieldSplit"
232 static PetscErrorCode PCSetUp_FieldSplit(PC pc)
233 {
234   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
235   PetscErrorCode    ierr;
236   PC_FieldSplitLink ilink;
237   PetscInt          i,nsplit,ccsize;
238   MatStructure      flag = pc->flag;
239   PetscTruth        sorted,getsub;
240 
241   PetscFunctionBegin;
242   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
243   nsplit = jac->nsplits;
244   ilink  = jac->head;
245 
246   /* get the matrices for each split */
247   if (!jac->issetup) {
248     PetscInt rstart,rend,nslots,bs;
249 
250     jac->issetup = PETSC_TRUE;
251 
252     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
253     bs     = jac->bs;
254     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
255     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
256     nslots = (rend - rstart)/bs;
257     for (i=0; i<nsplit; i++) {
258       if (jac->defaultsplit) {
259         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
260       } else if (!ilink->is) {
261         if (ilink->nfields > 1) {
262           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
263           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
264           for (j=0; j<nslots; j++) {
265             for (k=0; k<nfields; k++) {
266               ii[nfields*j + k] = rstart + bs*j + fields[k];
267             }
268           }
269           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
270           ierr = PetscFree(ii);CHKERRQ(ierr);
271         } else {
272           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
273         }
274       }
275       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
276       if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
277       ilink = ilink->next;
278     }
279   }
280 
281   ilink  = jac->head;
282   if (!jac->pmat) {
283     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
284     for (i=0; i<nsplit; i++) {
285       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
286       ilink = ilink->next;
287     }
288   } else {
289     for (i=0; i<nsplit; i++) {
290       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
291       ilink = ilink->next;
292     }
293   }
294   if (jac->realdiagonal) {
295     ilink = jac->head;
296     if (!jac->mat) {
297       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
298       for (i=0; i<nsplit; i++) {
299         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
300         ilink = ilink->next;
301       }
302     } else {
303       for (i=0; i<nsplit; i++) {
304         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
305         ilink = ilink->next;
306       }
307     }
308   } else {
309     jac->mat = jac->pmat;
310   }
311 
312   if (jac->type != PC_COMPOSITE_SCHUR) {
313     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
314     ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr);
315     if (getsub && jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
316       ilink  = jac->head;
317       if (!jac->Afield) {
318         ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
319         for (i=0; i<nsplit; i++) {
320           ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
321           ilink = ilink->next;
322         }
323       } else {
324         for (i=0; i<nsplit; i++) {
325           ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
326           ilink = ilink->next;
327         }
328       }
329     }
330   }
331 
332   if (jac->type == PC_COMPOSITE_SCHUR) {
333     IS       ccis;
334     PetscInt rstart,rend;
335     if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
336 
337     /* When extracting off-diagonal submatrices, we take complements from this range */
338     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
339 
340     /* need to handle case when one is resetting up the preconditioner */
341     if (jac->schur) {
342       ilink = jac->head;
343       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
344       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
345       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
346       ilink = ilink->next;
347       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
348       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
349       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
350       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
351       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
352 
353      } else {
354       KSP ksp;
355 
356       /* extract the B and C matrices */
357       ilink = jac->head;
358       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
359       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
360       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
361       ilink = ilink->next;
362       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
363       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
364       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
365       /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */
366       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
367       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
368       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
369       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
370 
371       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
372       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
373       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
374       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
375         PC pc;
376         ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
377         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
378         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
379       }
380       ierr  = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr);
381       ierr  = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr);
382       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
383       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
384 
385       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
386       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
387       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
388       ilink = jac->head;
389       ilink->x = jac->x[0]; ilink->y = jac->y[0];
390       ilink = ilink->next;
391       ilink->x = jac->x[1]; ilink->y = jac->y[1];
392     }
393   } else {
394     /* set up the individual PCs */
395     i    = 0;
396     ilink = jac->head;
397     while (ilink) {
398       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
399       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
400       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
401       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
402       i++;
403       ilink = ilink->next;
404     }
405 
406     /* create work vectors for each split */
407     if (!jac->x) {
408       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
409       ilink = jac->head;
410       for (i=0; i<nsplit; i++) {
411         Vec *vl,*vr;
412 
413         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
414         ilink->x  = *vr;
415         ilink->y  = *vl;
416         ierr      = PetscFree(vr);CHKERRQ(ierr);
417         ierr      = PetscFree(vl);CHKERRQ(ierr);
418         jac->x[i] = ilink->x;
419         jac->y[i] = ilink->y;
420         ilink     = ilink->next;
421       }
422     }
423   }
424 
425 
426   if (!jac->head->sctx) {
427     Vec xtmp;
428 
429     /* compute scatter contexts needed by multiplicative versions and non-default splits */
430 
431     ilink = jac->head;
432     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
433     for (i=0; i<nsplit; i++) {
434       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
435       ilink = ilink->next;
436     }
437     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
438   }
439   PetscFunctionReturn(0);
440 }
441 
442 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
443     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
444      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
445      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
446      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
447      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "PCApply_FieldSplit_Schur"
451 static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
452 {
453   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
454   PetscErrorCode    ierr;
455   KSP               ksp;
456   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
457 
458   PetscFunctionBegin;
459   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
460 
461   ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
462   ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
463   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
464   ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
465   ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
466   ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
467   ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
468 
469   ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
470   ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
471   ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
472 
473   ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
474   ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
475   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
476   ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
477   ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
478 
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "PCApply_FieldSplit"
484 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
485 {
486   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
487   PetscErrorCode    ierr;
488   PC_FieldSplitLink ilink = jac->head;
489   PetscInt          cnt;
490 
491   PetscFunctionBegin;
492   CHKMEMQ;
493   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
494   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
495 
496   if (jac->type == PC_COMPOSITE_ADDITIVE) {
497     if (jac->defaultsplit) {
498       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
499       while (ilink) {
500         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
501         ilink = ilink->next;
502       }
503       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
504     } else {
505       ierr = VecSet(y,0.0);CHKERRQ(ierr);
506       while (ilink) {
507         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
508         ilink = ilink->next;
509       }
510     }
511   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
512     if (!jac->w1) {
513       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
514       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
515     }
516     ierr = VecSet(y,0.0);CHKERRQ(ierr);
517     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
518     cnt = 1;
519     while (ilink->next) {
520       ilink = ilink->next;
521       if (jac->Afield) {
522         /* compute the residual only over the part of the vector needed */
523         ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
524         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
525         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
526         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
527         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
528         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
529         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530       } else {
531         /* compute the residual over the entire vector */
532 	ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
533 	ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
534 	ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
535       }
536     }
537     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
538       cnt -= 2;
539       while (ilink->previous) {
540         ilink = ilink->previous;
541         if (jac->Afield) {
542 	  /* compute the residual only over the part of the vector needed */
543 	  ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
544 	  ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
545 	  ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
546 	  ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
547 	  ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
548 	  ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
549 	  ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
550         } else {
551 	  ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
552 	  ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
553 	  ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
554         }
555       }
556     }
557   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
558   CHKMEMQ;
559   PetscFunctionReturn(0);
560 }
561 
562 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
563     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
564      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
565      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
566      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
567      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "PCApply_FieldSplit"
571 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
572 {
573   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
574   PetscErrorCode    ierr;
575   PC_FieldSplitLink ilink = jac->head;
576 
577   PetscFunctionBegin;
578   CHKMEMQ;
579   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
580   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
581 
582   if (jac->type == PC_COMPOSITE_ADDITIVE) {
583     if (jac->defaultsplit) {
584       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
585       while (ilink) {
586 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
587 	ilink = ilink->next;
588       }
589       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
590     } else {
591       ierr = VecSet(y,0.0);CHKERRQ(ierr);
592       while (ilink) {
593         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
594 	ilink = ilink->next;
595       }
596     }
597   } else {
598     if (!jac->w1) {
599       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
600       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
601     }
602     ierr = VecSet(y,0.0);CHKERRQ(ierr);
603     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
604       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
605       while (ilink->next) {
606         ilink = ilink->next;
607         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
608         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
609         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
610       }
611       while (ilink->previous) {
612         ilink = ilink->previous;
613         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
614         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
615         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
616       }
617     } else {
618       while (ilink->next) {   /* get to last entry in linked list */
619 	ilink = ilink->next;
620       }
621       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
622       while (ilink->previous) {
623 	ilink = ilink->previous;
624 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
625 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
626 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
627       }
628     }
629   }
630   CHKMEMQ;
631   PetscFunctionReturn(0);
632 }
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "PCDestroy_FieldSplit"
636 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
637 {
638   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
639   PetscErrorCode    ierr;
640   PC_FieldSplitLink ilink = jac->head,next;
641 
642   PetscFunctionBegin;
643   while (ilink) {
644     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
645     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
646     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
647     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
648     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
649     next = ilink->next;
650     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
651     ierr = PetscFree(ilink);CHKERRQ(ierr);
652     ilink = next;
653   }
654   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
655   if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);}
656   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
657   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
658   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
659   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
660   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
661   if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
662   if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);}
663   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
664   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
665   ierr = PetscFree(jac);CHKERRQ(ierr);
666   PetscFunctionReturn(0);
667 }
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
671 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
672 {
673   PetscErrorCode  ierr;
674   PetscInt        i = 0,nfields,*fields,bs;
675   PetscTruth      flg;
676   char            optionname[128];
677   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
678   PCCompositeType ctype;
679 
680   PetscFunctionBegin;
681   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
682   ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
683   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
684   if (flg) {
685     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
686   }
687 
688   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
689   if (flg) {
690     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
691   }
692 
693   /* Only setup fields once */
694   if ((jac->bs > 0) && (jac->nsplits == 0)) {
695     /* only allow user to set fields from command line if bs is already known.
696        otherwise user can set them in PCFieldSplitSetDefaults() */
697     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
698     while (PETSC_TRUE) {
699       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
700       nfields = jac->bs;
701       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
702       if (!flg) break;
703       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
704       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
705     }
706     ierr = PetscFree(fields);CHKERRQ(ierr);
707   }
708   ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,PETSC_NULL);CHKERRQ(ierr);
709   ierr = PetscOptionsTail();CHKERRQ(ierr);
710   PetscFunctionReturn(0);
711 }
712 
713 /*------------------------------------------------------------------------------------*/
714 
715 EXTERN_C_BEGIN
716 #undef __FUNCT__
717 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
718 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
719 {
720   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
721   PetscErrorCode    ierr;
722   PC_FieldSplitLink ilink,next = jac->head;
723   char              prefix[128];
724   PetscInt          i;
725 
726   PetscFunctionBegin;
727   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
728   for (i=0; i<n; i++) {
729     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
730     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
731   }
732   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
733   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
734   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
735   ilink->nfields = n;
736   ilink->next    = PETSC_NULL;
737   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
738   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
739   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
740 
741   if (((PetscObject)pc)->prefix) {
742     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
743   } else {
744     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
745   }
746   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
747 
748   if (!next) {
749     jac->head       = ilink;
750     ilink->previous = PETSC_NULL;
751   } else {
752     while (next->next) {
753       next = next->next;
754     }
755     next->next      = ilink;
756     ilink->previous = next;
757   }
758   jac->nsplits++;
759   PetscFunctionReturn(0);
760 }
761 EXTERN_C_END
762 
763 EXTERN_C_BEGIN
764 #undef __FUNCT__
765 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
766 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
767 {
768   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
769   PetscErrorCode ierr;
770 
771   PetscFunctionBegin;
772   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
773   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
774   (*subksp)[1] = jac->kspschur;
775   *n = jac->nsplits;
776   PetscFunctionReturn(0);
777 }
778 EXTERN_C_END
779 
780 EXTERN_C_BEGIN
781 #undef __FUNCT__
782 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
783 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
784 {
785   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
786   PetscErrorCode    ierr;
787   PetscInt          cnt = 0;
788   PC_FieldSplitLink ilink = jac->head;
789 
790   PetscFunctionBegin;
791   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
792   while (ilink) {
793     (*subksp)[cnt++] = ilink->ksp;
794     ilink = ilink->next;
795   }
796   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
797   *n = jac->nsplits;
798   PetscFunctionReturn(0);
799 }
800 EXTERN_C_END
801 
802 EXTERN_C_BEGIN
803 #undef __FUNCT__
804 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
805 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
806 {
807   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
808   PetscErrorCode    ierr;
809   PC_FieldSplitLink ilink, next = jac->head;
810   char              prefix[128];
811 
812   PetscFunctionBegin;
813   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
814   ilink->is      = is;
815   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
816   ilink->next    = PETSC_NULL;
817   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
818   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
819   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
820 
821   if (((PetscObject)pc)->prefix) {
822     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
823   } else {
824     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
825   }
826   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
827 
828   if (!next) {
829     jac->head       = ilink;
830     ilink->previous = PETSC_NULL;
831   } else {
832     while (next->next) {
833       next = next->next;
834     }
835     next->next      = ilink;
836     ilink->previous = next;
837   }
838   jac->nsplits++;
839 
840   PetscFunctionReturn(0);
841 }
842 EXTERN_C_END
843 
844 #undef __FUNCT__
845 #define __FUNCT__ "PCFieldSplitSetFields"
846 /*@
847     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
848 
849     Collective on PC
850 
851     Input Parameters:
852 +   pc  - the preconditioner context
853 .   n - the number of fields in this split
854 .   fields - the fields in this split
855 
856     Level: intermediate
857 
858     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
859 
860      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
861      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
862      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
863      where the numbered entries indicate what is in the field.
864 
865 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
866 
867 @*/
868 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
869 {
870   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
874   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
875   if (f) {
876     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
877   }
878   PetscFunctionReturn(0);
879 }
880 
881 #undef __FUNCT__
882 #define __FUNCT__ "PCFieldSplitSetIS"
883 /*@
884     PCFieldSplitSetIS - Sets the exact elements for field
885 
886     Collective on PC
887 
888     Input Parameters:
889 +   pc  - the preconditioner context
890 .   is - the index set that defines the vector elements in this field
891 
892 
893     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
894 
895     Level: intermediate
896 
897 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
898 
899 @*/
900 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
901 {
902   PetscErrorCode ierr,(*f)(PC,IS);
903 
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
906   PetscValidHeaderSpecific(is,IS_COOKIE,1);
907   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
908   if (f) {
909     ierr = (*f)(pc,is);CHKERRQ(ierr);
910   }
911   PetscFunctionReturn(0);
912 }
913 
914 #undef __FUNCT__
915 #define __FUNCT__ "PCFieldSplitSetBlockSize"
916 /*@
917     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
918       fieldsplit preconditioner. If not set the matrix block size is used.
919 
920     Collective on PC
921 
922     Input Parameters:
923 +   pc  - the preconditioner context
924 -   bs - the block size
925 
926     Level: intermediate
927 
928 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
929 
930 @*/
931 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
932 {
933   PetscErrorCode ierr,(*f)(PC,PetscInt);
934 
935   PetscFunctionBegin;
936   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
937   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
938   if (f) {
939     ierr = (*f)(pc,bs);CHKERRQ(ierr);
940   }
941   PetscFunctionReturn(0);
942 }
943 
944 #undef __FUNCT__
945 #define __FUNCT__ "PCFieldSplitGetSubKSP"
946 /*@C
947    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
948 
949    Collective on KSP
950 
951    Input Parameter:
952 .  pc - the preconditioner context
953 
954    Output Parameters:
955 +  n - the number of split
956 -  pc - the array of KSP contexts
957 
958    Note:
959    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
960    (not the KSP just the array that contains them).
961 
962    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
963 
964    Level: advanced
965 
966 .seealso: PCFIELDSPLIT
967 @*/
968 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
969 {
970   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
971 
972   PetscFunctionBegin;
973   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
974   PetscValidIntPointer(n,2);
975   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
976   if (f) {
977     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
978   } else {
979     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
980   }
981   PetscFunctionReturn(0);
982 }
983 
984 #undef __FUNCT__
985 #define __FUNCT__ "PCFieldSplitSchurPrecondition"
986 /*@
987     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
988       D matrix. Otherwise no preconditioner is used.
989 
990     Collective on PC
991 
992     Input Parameters:
993 +   pc  - the preconditioner context
994 .   ptype - which matrix to use for preconditioning the Schur complement
995 -   userpre - matrix to use for preconditioning, or PETSC_NULL
996 
997     Notes:
998     The default is to use the block on the diagonal of the preconditioning matrix.  This is D, in the (1,1) position.
999     There are currently no preconditioners that work directly with the Schur complement so setting
1000     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
1001 
1002     Options Database:
1003 .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1004 
1005     Level: intermediate
1006 
1007 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1008 
1009 @*/
1010 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1011 {
1012   PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat);
1013 
1014   PetscFunctionBegin;
1015   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1016   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr);
1017   if (f) {
1018     ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr);
1019   }
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 EXTERN_C_BEGIN
1024 #undef __FUNCT__
1025 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1026 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1027 {
1028   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1029   PetscErrorCode  ierr;
1030 
1031   PetscFunctionBegin;
1032   jac->schurpre = ptype;
1033   if (pre) {
1034     if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
1035     jac->schur_user = pre;
1036     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1037   }
1038   PetscFunctionReturn(0);
1039 }
1040 EXTERN_C_END
1041 
1042 #undef __FUNCT__
1043 #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
1044 /*@C
1045    PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement
1046 
1047    Collective on KSP
1048 
1049    Input Parameter:
1050 .  pc - the preconditioner context
1051 
1052    Output Parameters:
1053 +  A - the (0,0) block
1054 .  B - the (0,1) block
1055 .  C - the (1,0) block
1056 -  D - the (1,1) block
1057 
1058    Level: advanced
1059 
1060 .seealso: PCFIELDSPLIT
1061 @*/
1062 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D)
1063 {
1064   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
1065 
1066   PetscFunctionBegin;
1067   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1068   if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");}
1069   if (A) *A = jac->pmat[0];
1070   if (B) *B = jac->B;
1071   if (C) *C = jac->C;
1072   if (D) *D = jac->pmat[1];
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 EXTERN_C_BEGIN
1077 #undef __FUNCT__
1078 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1079 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
1080 {
1081   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1082   PetscErrorCode ierr;
1083 
1084   PetscFunctionBegin;
1085   jac->type = type;
1086   if (type == PC_COMPOSITE_SCHUR) {
1087     pc->ops->apply = PCApply_FieldSplit_Schur;
1088     pc->ops->view  = PCView_FieldSplit_Schur;
1089     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1090     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1091 
1092   } else {
1093     pc->ops->apply = PCApply_FieldSplit;
1094     pc->ops->view  = PCView_FieldSplit;
1095     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1096     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1097   }
1098   PetscFunctionReturn(0);
1099 }
1100 EXTERN_C_END
1101 
1102 EXTERN_C_BEGIN
1103 #undef __FUNCT__
1104 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
1105 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
1106 {
1107   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1108 
1109   PetscFunctionBegin;
1110   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1111   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
1112   jac->bs = bs;
1113   PetscFunctionReturn(0);
1114 }
1115 EXTERN_C_END
1116 
1117 #undef __FUNCT__
1118 #define __FUNCT__ "PCFieldSplitSetType"
1119 /*@
1120    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
1121 
1122    Collective on PC
1123 
1124    Input Parameter:
1125 .  pc - the preconditioner context
1126 .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1127 
1128    Options Database Key:
1129 .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
1130 
1131    Level: Developer
1132 
1133 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1134 
1135 .seealso: PCCompositeSetType()
1136 
1137 @*/
1138 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
1139 {
1140   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
1141 
1142   PetscFunctionBegin;
1143   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1144   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
1145   if (f) {
1146     ierr = (*f)(pc,type);CHKERRQ(ierr);
1147   }
1148   PetscFunctionReturn(0);
1149 }
1150 
1151 /* -------------------------------------------------------------------------------------*/
1152 /*MC
1153    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1154                   fields or groups of fields
1155 
1156 
1157      To set options on the solvers for each block append -fieldsplit_ to all the PC
1158         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
1159 
1160      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
1161          and set the options directly on the resulting KSP object
1162 
1163    Level: intermediate
1164 
1165    Options Database Keys:
1166 +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
1167 .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
1168                               been supplied explicitly by -pc_fieldsplit_%d_fields
1169 .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
1170 .    -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1171 .    -pc_fieldsplit_schur_precondition <true,false> default is true
1172 
1173 -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
1174      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
1175 
1176 
1177    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1178      to define a field by an arbitrary collection of entries.
1179 
1180       If no fields are set the default is used. The fields are defined by entries strided by bs,
1181       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1182       if this is not called the block size defaults to the blocksize of the second matrix passed
1183       to KSPSetOperators()/PCSetOperators().
1184 
1185       Currently for the multiplicative version, the updated residual needed for the next field
1186      solve is computed via a matrix vector product over the entire array. An optimization would be
1187      to update the residual only for the part of the right hand side associated with the next field
1188      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1189      part of the matrix needed to just update part of the residual).
1190 
1191      For the Schur complement preconditioner if J = ( A B )
1192                                                     ( C D )
1193      the preconditioner is
1194               (I   -B inv(A)) ( inv(A)   0    ) (I         0  )
1195               (0    I       ) (   0    inv(S) ) (-C inv(A) I  )
1196      where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1197      inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is
1198      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1199      D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1200      option.
1201 
1202      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1203      is used automatically for a second block.
1204 
1205    Concepts: physics based preconditioners, block preconditioners
1206 
1207 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners
1208            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1209 M*/
1210 
1211 EXTERN_C_BEGIN
1212 #undef __FUNCT__
1213 #define __FUNCT__ "PCCreate_FieldSplit"
1214 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
1215 {
1216   PetscErrorCode ierr;
1217   PC_FieldSplit  *jac;
1218 
1219   PetscFunctionBegin;
1220   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
1221   jac->bs        = -1;
1222   jac->nsplits   = 0;
1223   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1224   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1225 
1226   pc->data     = (void*)jac;
1227 
1228   pc->ops->apply             = PCApply_FieldSplit;
1229   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
1230   pc->ops->setup             = PCSetUp_FieldSplit;
1231   pc->ops->destroy           = PCDestroy_FieldSplit;
1232   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
1233   pc->ops->view              = PCView_FieldSplit;
1234   pc->ops->applyrichardson   = 0;
1235 
1236   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
1237                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1238   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
1239                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1240   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1241                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
1242   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
1243                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
1244   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
1245                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
1246   PetscFunctionReturn(0);
1247 }
1248 EXTERN_C_END
1249 
1250 
1251 /*MC
1252    Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an
1253           overview of these methods.
1254 
1255       Consider the solution to ( A B ) (x_1)  =  (b_1)
1256                                ( C D ) (x_2)     (b_2)
1257 
1258       Important special cases, the Stokes equation: C = B' and D = 0  (A   B) (x_1) = (b_1)
1259                                                                        B'  0) (x_2)   (b_2)
1260 
1261       One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners
1262       for this block system.
1263 
1264       Consider an additional matrix (Ap  Bp)
1265                                     (Cp  Dp) where some or all of the entries may be the same as
1266       in the original matrix (for example Ap == A).
1267 
1268       In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the
1269       approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...)
1270 
1271       Block Jacobi:   x_1 = A^ b_1
1272                       x_2 = D^ b_2
1273 
1274       Lower block Gauss-Seidel:   x_1 = A^ b_1
1275                             x_2 = D^ (b_2 - C x_1)       variant x_2 = D^ (b_2 - Cp x_1)
1276 
1277       Symmetric Gauss-Seidel:  x_1 = x_1 + A^(b_1 - A x_1 - B x_2)    variant  x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2)
1278           Interestingly this form is not actually a symmetric matrix, the symmetric version is
1279                               x_1 = A^(b_1 - B x_2)      variant x_1 = A^(b_1 - Bp x_2)
1280 
1281    Level: intermediate
1282 
1283 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT
1284            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1285 M*/
1286