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