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