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