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