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