xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 4d0a8057dc766494e73f985a8c029a3f3edde3ad)
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   PC_FieldSplitLink head;
30 } PC_FieldSplit;
31 
32 /*
33     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
34    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
35    PC you could change this.
36 */
37 #undef __FUNCT__
38 #define __FUNCT__ "PCView_FieldSplit"
39 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
40 {
41   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
42   PetscErrorCode    ierr;
43   PetscTruth        iascii;
44   PetscInt          i,j;
45   PC_FieldSplitLink ilink = jac->head;
46 
47   PetscFunctionBegin;
48   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
49   if (iascii) {
50     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
51     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
52     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
53     for (i=0; i<jac->nsplits; i++) {
54       if (ilink->fields) {
55 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
56 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
57 	for (j=0; j<ilink->nfields; j++) {
58 	  if (j > 0) {
59 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
60 	  }
61 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
62 	}
63 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
64         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
65       } else {
66 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
67       }
68       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
69       ilink = ilink->next;
70     }
71     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
72   } else {
73     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
74   }
75   PetscFunctionReturn(0);
76 }
77 
78 #undef __FUNCT__
79 #define __FUNCT__ "PCFieldSplitSetDefaults"
80 static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
81 {
82   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
83   PetscErrorCode    ierr;
84   PC_FieldSplitLink ilink = jac->head;
85   PetscInt          i = 0,*ifields,nfields;
86   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
87   char              optionname[128];
88 
89   PetscFunctionBegin;
90   if (!ilink) {
91 
92     if (jac->bs <= 0) {
93       if (pc->pmat) {
94         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
95       } else {
96         jac->bs = 1;
97       }
98     }
99 
100     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
101     if (!flg) {
102       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
103          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
104       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
105       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
106       while (PETSC_TRUE) {
107         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
108         nfields = jac->bs;
109         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
110         if (!flg2) break;
111         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
112         flg = PETSC_FALSE;
113         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
114       }
115       ierr = PetscFree(ifields);CHKERRQ(ierr);
116     }
117 
118     if (flg) {
119       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
120       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
121       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
122       while (ilink) {
123 	for (i=0; i<ilink->nfields; i++) {
124 	  fields[ilink->fields[i]] = PETSC_TRUE;
125 	}
126 	ilink = ilink->next;
127       }
128       jac->defaultsplit = PETSC_TRUE;
129       for (i=0; i<jac->bs; i++) {
130 	if (!fields[i]) {
131 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
132 	} else {
133 	  jac->defaultsplit = PETSC_FALSE;
134 	}
135       }
136       ierr = PetscFree(fields);CHKERRQ(ierr);
137     }
138   }
139   PetscFunctionReturn(0);
140 }
141 
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "PCSetUp_FieldSplit"
145 static PetscErrorCode PCSetUp_FieldSplit(PC pc)
146 {
147   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
148   PetscErrorCode    ierr;
149   PC_FieldSplitLink ilink;
150   PetscInt          i,nsplit,ccsize;
151   MatStructure      flag = pc->flag;
152   PetscTruth        sorted,getsub;
153 
154   PetscFunctionBegin;
155   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
156   nsplit = jac->nsplits;
157   ilink  = jac->head;
158 
159   /* get the matrices for each split */
160   if (!jac->issetup) {
161     PetscInt rstart,rend,nslots,bs;
162 
163     jac->issetup = PETSC_TRUE;
164 
165     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
166     bs     = jac->bs;
167     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
168     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
169     nslots = (rend - rstart)/bs;
170     for (i=0; i<nsplit; i++) {
171       if (jac->defaultsplit) {
172         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
173         ilink->csize = ccsize/nsplit;
174       } else if (!ilink->is) {
175         if (ilink->nfields > 1) {
176           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
177           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
178           for (j=0; j<nslots; j++) {
179             for (k=0; k<nfields; k++) {
180               ii[nfields*j + k] = rstart + bs*j + fields[k];
181             }
182           }
183           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
184           ierr = PetscFree(ii);CHKERRQ(ierr);
185         } else {
186           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
187         }
188         ilink->csize = (ccsize/bs)*ilink->nfields;
189         ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
190         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
191       }
192       ierr  = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
193       ilink = ilink->next;
194     }
195   }
196 
197   ilink  = jac->head;
198   if (!jac->pmat) {
199     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
200     for (i=0; i<nsplit; i++) {
201       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
202       ilink = ilink->next;
203     }
204   } else {
205     for (i=0; i<nsplit; i++) {
206       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
207       ilink = ilink->next;
208     }
209   }
210 
211   /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
212   ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr);
213   if (getsub && jac->type != PC_COMPOSITE_ADDITIVE) {
214     ilink  = jac->head;
215     if (!jac->Afield) {
216       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
217       for (i=0; i<nsplit; i++) {
218 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,ilink->csize,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
219 	ilink = ilink->next;
220       }
221     } else {
222       for (i=0; i<nsplit; i++) {
223 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,ilink->csize,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
224 	ilink = ilink->next;
225       }
226     }
227   }
228 
229 
230   /* set up the individual PCs */
231   i    = 0;
232   ilink = jac->head;
233   while (ilink) {
234     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
235     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
236     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
237     i++;
238     ilink = ilink->next;
239   }
240 
241   /* create work vectors for each split */
242   if (!jac->x) {
243     Vec xtmp;
244     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
245     ilink = jac->head;
246     for (i=0; i<nsplit; i++) {
247       Vec *vl,*vr;
248 
249       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
250       ilink->x  = *vr;
251       ilink->y  = *vl;
252       ierr      = PetscFree(vr);CHKERRQ(ierr);
253       ierr      = PetscFree(vl);CHKERRQ(ierr);
254       jac->x[i] = ilink->x;
255       jac->y[i] = ilink->y;
256       ilink     = ilink->next;
257     }
258     /* compute scatter contexts needed by multiplicative versions and non-default splits */
259 
260     ilink = jac->head;
261     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
262     for (i=0; i<nsplit; i++) {
263       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
264       ilink = ilink->next;
265     }
266     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
267   }
268   PetscFunctionReturn(0);
269 }
270 
271 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
272     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
273      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
274      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
275      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
276      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "PCApply_FieldSplit"
280 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
281 {
282   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
283   PetscErrorCode    ierr;
284   PC_FieldSplitLink ilink = jac->head;
285   PetscInt          bs;
286 
287   PetscFunctionBegin;
288   CHKMEMQ;
289   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
290   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
291   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
292 
293   if (jac->type == PC_COMPOSITE_ADDITIVE) {
294     if (jac->defaultsplit) {
295       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
296       while (ilink) {
297         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
298         ilink = ilink->next;
299       }
300       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
301     } else {
302       ierr = VecSet(y,0.0);CHKERRQ(ierr);
303       while (ilink) {
304         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
305         ilink = ilink->next;
306       }
307     }
308   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
309     if (!jac->w1) {
310       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
311       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
312     }
313     ierr = VecSet(y,0.0);CHKERRQ(ierr);
314     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
315     while (ilink->next) {
316       ilink = ilink->next;
317       ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
318       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
319       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
320     }
321     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
322       while (ilink->previous) {
323         ilink = ilink->previous;
324         ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
325         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
326         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
327       }
328     }
329   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
330   CHKMEMQ;
331   PetscFunctionReturn(0);
332 }
333 
334 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
335     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
336      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
337      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
338      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
339      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "PCApply_FieldSplit"
343 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
344 {
345   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
346   PetscErrorCode    ierr;
347   PC_FieldSplitLink ilink = jac->head;
348   PetscInt          bs;
349 
350   PetscFunctionBegin;
351   CHKMEMQ;
352   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
353   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
354   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
355 
356   if (jac->type == PC_COMPOSITE_ADDITIVE) {
357     if (jac->defaultsplit) {
358       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
359       while (ilink) {
360 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
361 	ilink = ilink->next;
362       }
363       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
364     } else {
365       ierr = VecSet(y,0.0);CHKERRQ(ierr);
366       while (ilink) {
367         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
368 	ilink = ilink->next;
369       }
370     }
371   } else {
372     if (!jac->w1) {
373       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
374       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
375     }
376     ierr = VecSet(y,0.0);CHKERRQ(ierr);
377     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
378       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
379       while (ilink->next) {
380         ilink = ilink->next;
381         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
382         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
383         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
384       }
385       while (ilink->previous) {
386         ilink = ilink->previous;
387         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
388         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
389         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
390       }
391     } else {
392       while (ilink->next) {   /* get to last entry in linked list */
393 	ilink = ilink->next;
394       }
395       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
396       while (ilink->previous) {
397 	ilink = ilink->previous;
398 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
399 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
400 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
401       }
402     }
403   }
404   CHKMEMQ;
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "PCDestroy_FieldSplit"
410 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
411 {
412   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
413   PetscErrorCode    ierr;
414   PC_FieldSplitLink ilink = jac->head,next;
415 
416   PetscFunctionBegin;
417   while (ilink) {
418     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
419     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
420     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
421     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
422     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
423     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
424     next = ilink->next;
425     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
426     ierr = PetscFree(ilink);CHKERRQ(ierr);
427     ilink = next;
428   }
429   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
430   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
431   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
432   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
433   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
434   ierr = PetscFree(jac);CHKERRQ(ierr);
435   PetscFunctionReturn(0);
436 }
437 
438 #undef __FUNCT__
439 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
440 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
441 {
442   PetscErrorCode ierr;
443   PetscInt       i = 0,nfields,*fields,bs;
444   PetscTruth     flg;
445   char           optionname[128];
446   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
447 
448   PetscFunctionBegin;
449   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
450   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
451   if (flg) {
452     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
453   }
454 
455   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
456 
457   if (jac->bs > 0) {
458     /* only allow user to set fields from command line if bs is already known.
459        otherwise user can set them in PCFieldSplitSetDefaults() */
460     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
461     while (PETSC_TRUE) {
462       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
463       nfields = jac->bs;
464       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
465       if (!flg) break;
466       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
467       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
468     }
469     ierr = PetscFree(fields);CHKERRQ(ierr);
470   }
471   ierr = PetscOptionsTail();CHKERRQ(ierr);
472   PetscFunctionReturn(0);
473 }
474 
475 /*------------------------------------------------------------------------------------*/
476 
477 EXTERN_C_BEGIN
478 #undef __FUNCT__
479 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
480 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
481 {
482   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
483   PetscErrorCode    ierr;
484   PC_FieldSplitLink ilink,next = jac->head;
485   char              prefix[128];
486   PetscInt          i;
487 
488   PetscFunctionBegin;
489   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
490   for (i=0; i<n; i++) {
491     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
492     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
493   }
494   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
495   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
496   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
497   ilink->nfields = n;
498   ilink->next    = PETSC_NULL;
499   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
500   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
501 
502   if (((PetscObject)pc)->prefix) {
503     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
504   } else {
505     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
506   }
507   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
508 
509   if (!next) {
510     jac->head       = ilink;
511     ilink->previous = PETSC_NULL;
512   } else {
513     while (next->next) {
514       next = next->next;
515     }
516     next->next      = ilink;
517     ilink->previous = next;
518   }
519   jac->nsplits++;
520   PetscFunctionReturn(0);
521 }
522 EXTERN_C_END
523 
524 
525 EXTERN_C_BEGIN
526 #undef __FUNCT__
527 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
528 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
529 {
530   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
531   PetscErrorCode    ierr;
532   PetscInt          cnt = 0;
533   PC_FieldSplitLink ilink = jac->head;
534 
535   PetscFunctionBegin;
536   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
537   while (ilink) {
538     (*subksp)[cnt++] = ilink->ksp;
539     ilink = ilink->next;
540   }
541   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
542   *n = jac->nsplits;
543   PetscFunctionReturn(0);
544 }
545 EXTERN_C_END
546 
547 EXTERN_C_BEGIN
548 #undef __FUNCT__
549 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
550 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
551 {
552   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
553   PetscErrorCode    ierr;
554     PC_FieldSplitLink ilink, next = jac->head;
555   char              prefix[128];
556 
557   PetscFunctionBegin;
558   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
559   ilink->is      = is;
560   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
561   ilink->next    = PETSC_NULL;
562   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
563   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
564 
565   if (((PetscObject)pc)->prefix) {
566     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
567   } else {
568     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
569   }
570   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
571 
572   if (!next) {
573     jac->head       = ilink;
574     ilink->previous = PETSC_NULL;
575   } else {
576     while (next->next) {
577       next = next->next;
578     }
579     next->next      = ilink;
580     ilink->previous = next;
581   }
582   jac->nsplits++;
583 
584   PetscFunctionReturn(0);
585 }
586 EXTERN_C_END
587 
588 #undef __FUNCT__
589 #define __FUNCT__ "PCFieldSplitSetFields"
590 /*@
591     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
592 
593     Collective on PC
594 
595     Input Parameters:
596 +   pc  - the preconditioner context
597 .   n - the number of fields in this split
598 .   fields - the fields in this split
599 
600     Level: intermediate
601 
602     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
603 
604      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
605      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
606      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
607      where the numbered entries indicate what is in the field.
608 
609 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
610 
611 @*/
612 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
613 {
614   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
615 
616   PetscFunctionBegin;
617   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
618   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
619   if (f) {
620     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
621   }
622   PetscFunctionReturn(0);
623 }
624 
625 #undef __FUNCT__
626 #define __FUNCT__ "PCFieldSplitSetIS"
627 /*@
628     PCFieldSplitSetIS - Sets the exact elements for field
629 
630     Collective on PC
631 
632     Input Parameters:
633 +   pc  - the preconditioner context
634 .   is - the index set that defines the vector elements in this field
635 
636 
637     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
638 
639     Level: intermediate
640 
641 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
642 
643 @*/
644 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
645 {
646   PetscErrorCode ierr,(*f)(PC,IS);
647 
648   PetscFunctionBegin;
649   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
650   PetscValidHeaderSpecific(is,IS_COOKIE,1);
651   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
652   if (f) {
653     ierr = (*f)(pc,is);CHKERRQ(ierr);
654   }
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "PCFieldSplitSetBlockSize"
660 /*@
661     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
662       fieldsplit preconditioner. If not set the matrix block size is used.
663 
664     Collective on PC
665 
666     Input Parameters:
667 +   pc  - the preconditioner context
668 -   bs - the block size
669 
670     Level: intermediate
671 
672 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
673 
674 @*/
675 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
676 {
677   PetscErrorCode ierr,(*f)(PC,PetscInt);
678 
679   PetscFunctionBegin;
680   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
681   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
682   if (f) {
683     ierr = (*f)(pc,bs);CHKERRQ(ierr);
684   }
685   PetscFunctionReturn(0);
686 }
687 
688 #undef __FUNCT__
689 #define __FUNCT__ "PCFieldSplitGetSubKSP"
690 /*@C
691    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
692 
693    Collective on KSP
694 
695    Input Parameter:
696 .  pc - the preconditioner context
697 
698    Output Parameters:
699 +  n - the number of split
700 -  pc - the array of KSP contexts
701 
702    Note:
703    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
704    (not the KSP just the array that contains them).
705 
706    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
707 
708    Level: advanced
709 
710 .seealso: PCFIELDSPLIT
711 @*/
712 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
713 {
714   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
715 
716   PetscFunctionBegin;
717   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
718   PetscValidIntPointer(n,2);
719   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
720   if (f) {
721     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
722   } else {
723     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
724   }
725   PetscFunctionReturn(0);
726 }
727 
728 EXTERN_C_BEGIN
729 #undef __FUNCT__
730 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
731 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
732 {
733   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
734 
735   PetscFunctionBegin;
736   jac->type = type;
737   PetscFunctionReturn(0);
738 }
739 EXTERN_C_END
740 
741 EXTERN_C_BEGIN
742 #undef __FUNCT__
743 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
744 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
745 {
746   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
747 
748   PetscFunctionBegin;
749   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
750   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);
751   jac->bs = bs;
752   PetscFunctionReturn(0);
753 }
754 EXTERN_C_END
755 
756 #undef __FUNCT__
757 #define __FUNCT__ "PCFieldSplitSetType"
758 /*@
759    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
760 
761    Collective on PC
762 
763    Input Parameter:
764 .  pc - the preconditioner context
765 .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
766 
767    Options Database Key:
768 .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
769 
770    Level: Developer
771 
772 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
773 
774 .seealso: PCCompositeSetType()
775 
776 @*/
777 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
778 {
779   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
780 
781   PetscFunctionBegin;
782   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
783   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
784   if (f) {
785     ierr = (*f)(pc,type);CHKERRQ(ierr);
786   }
787   PetscFunctionReturn(0);
788 }
789 
790 /* -------------------------------------------------------------------------------------*/
791 /*MC
792    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
793                   fields or groups of fields
794 
795 
796      To set options on the solvers for each block append -sub_ to all the PC
797         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
798 
799      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
800          and set the options directly on the resulting KSP object
801 
802    Level: intermediate
803 
804    Options Database Keys:
805 +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
806 .   -pc_splitfield_default - automatically add any fields to additional splits that have not
807                               been supplied explicitly by -pc_splitfield_%d_fields
808 .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
809 -   -pc_splitfield_type <additive,multiplicative>
810 
811    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
812      to define a field by an arbitrary collection of entries.
813 
814       If no fields are set the default is used. The fields are defined by entries strided by bs,
815       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
816       if this is not called the block size defaults to the blocksize of the second matrix passed
817       to KSPSetOperators()/PCSetOperators().
818 
819       Currently for the multiplicative version, the updated residual needed for the next field
820      solve is computed via a matrix vector product over the entire array. An optimization would be
821      to update the residual only for the part of the right hand side associated with the next field
822      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
823      part of the matrix needed to just update part of the residual).
824 
825    Concepts: physics based preconditioners
826 
827 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
828            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS()
829 M*/
830 
831 EXTERN_C_BEGIN
832 #undef __FUNCT__
833 #define __FUNCT__ "PCCreate_FieldSplit"
834 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
835 {
836   PetscErrorCode ierr;
837   PC_FieldSplit  *jac;
838 
839   PetscFunctionBegin;
840   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
841   jac->bs        = -1;
842   jac->nsplits   = 0;
843   jac->type      = PC_COMPOSITE_ADDITIVE;
844 
845   pc->data     = (void*)jac;
846 
847   pc->ops->apply             = PCApply_FieldSplit;
848   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
849   pc->ops->setup             = PCSetUp_FieldSplit;
850   pc->ops->destroy           = PCDestroy_FieldSplit;
851   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
852   pc->ops->view              = PCView_FieldSplit;
853   pc->ops->applyrichardson   = 0;
854 
855   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
856                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
857   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
858                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
859   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
860                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
861   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
862                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
863   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
864                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
865   PetscFunctionReturn(0);
866 }
867 EXTERN_C_END
868 
869 
870