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