xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision ccb205f84456404d0833d80f1cbd41e068105bb1)
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   PC_FieldSplitLink next,previous;
16 };
17 
18 typedef struct {
19   PCCompositeType   type;              /* additive or multiplicative */
20   PetscTruth        defaultsplit;
21   PetscInt          bs;
22   PetscInt          nsplits,*csize;
23   Vec               *x,*y,w1,w2;
24   Mat               *pmat;
25   IS                *is,*cis;
26   PC_FieldSplitLink head;
27 } PC_FieldSplit;
28 
29 #undef __FUNCT__
30 #define __FUNCT__ "PCView_FieldSplit"
31 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
32 {
33   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
34   PetscErrorCode    ierr;
35   PetscTruth        iascii;
36   PetscInt          i,j;
37   PC_FieldSplitLink ilink = jac->head;
38 
39   PetscFunctionBegin;
40   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
41   if (iascii) {
42     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
43     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
44     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
45     for (i=0; i<jac->nsplits; i++) {
46       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
47       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
48       for (j=0; j<ilink->nfields; j++) {
49         if (j > 0) {
50           ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
51         }
52         ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
53       }
54       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
55       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
56       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
57       ilink = ilink->next;
58     }
59     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
60   } else {
61     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
62   }
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "PCFieldSplitSetDefaults"
68 static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
69 {
70   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
71   PetscErrorCode    ierr;
72   PC_FieldSplitLink ilink = jac->head;
73   PetscInt          i;
74   PetscTruth        flg = PETSC_FALSE,*fields;
75 
76   PetscFunctionBegin;
77   ierr = PetscOptionsGetTruth(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
78   if (!ilink || flg) {
79     ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
80     if (jac->bs <= 0) {
81       ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
82     }
83     ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
84     ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
85     while (ilink) {
86       for (i=0; i<ilink->nfields; i++) {
87         fields[ilink->fields[i]] = PETSC_TRUE;
88       }
89       ilink = ilink->next;
90     }
91     jac->defaultsplit = PETSC_TRUE;
92     for (i=0; i<jac->bs; i++) {
93       if (!fields[i]) {
94 	ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
95       } else {
96         jac->defaultsplit = PETSC_FALSE;
97       }
98     }
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "PCSetUp_FieldSplit"
106 static PetscErrorCode PCSetUp_FieldSplit(PC pc)
107 {
108   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
109   PetscErrorCode    ierr;
110   PC_FieldSplitLink ilink;
111   PetscInt          i,nsplit,ccsize;
112   MatStructure      flag = pc->flag;
113   PetscTruth        sorted;
114 
115   PetscFunctionBegin;
116   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
117   nsplit = jac->nsplits;
118   ilink  = jac->head;
119 
120   /* get the matrices for each split */
121   if (!jac->is) {
122     PetscInt rstart,rend,nslots,bs;
123 
124     bs     = jac->bs;
125     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
126     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
127     nslots = (rend - rstart)/bs;
128     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr);
129     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->cis);CHKERRQ(ierr);
130     ierr   = PetscMalloc(nsplit*sizeof(PetscInt),&jac->csize);CHKERRQ(ierr);
131     for (i=0; i<nsplit; i++) {
132       if (jac->defaultsplit) {
133 	ierr     = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr);
134         jac->csize[i] = ccsize/nsplit;
135       } else {
136 	printf("ilink nfields %d\n",ilink->nfields);
137         if (ilink->nfields > 1) {
138 	  PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
139 	  ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
140 	  for (j=0; j<nslots; j++) {
141 	    for (k=0; k<nfields; k++) {
142 	      ii[nfields*j + k] = rstart + bs*j + fields[k];
143 	    }
144 	  }
145 	  ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr);
146 	  ierr = PetscFree(ii);CHKERRQ(ierr);
147         } else {
148           ierr = ISCreateStride(pc->comm,nslots,ilink->fields[0],bs,&jac->is[i]);CHKERRQ(ierr);
149         }
150         jac->csize[i] = (ccsize/bs)*ilink->nfields;
151         ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr);
152         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
153         ilink = ilink->next;
154       }
155       ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr);
156     }
157   }
158 
159   if (!jac->pmat) {
160     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
161     for (i=0; i<nsplit; i++) {
162       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
163     }
164   } else {
165     for (i=0; i<nsplit; i++) {
166       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
167     }
168   }
169 
170   /* set up the individual PCs */
171   i    = 0;
172   ilink = jac->head;
173   while (ilink) {
174     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
175     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
176     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
177     i++;
178     ilink = ilink->next;
179   }
180 
181   /* create work vectors for each split */
182   if (!jac->x) {
183     Vec xtmp;
184     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
185     ilink = jac->head;
186     for (i=0; i<nsplit; i++) {
187       Vec *vl,*vr;
188 
189       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
190       ilink->x  = *vr;
191       ilink->y  = *vl;
192       ierr      = PetscFree(vr);CHKERRQ(ierr);
193       ierr      = PetscFree(vl);CHKERRQ(ierr);
194       jac->x[i] = ilink->x;
195       jac->y[i] = ilink->y;
196       ilink     = ilink->next;
197     }
198     /* compute scatter contexts needed by multiplicative versions and non-default splits */
199 
200     ilink = jac->head;
201     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
202     for (i=0; i<nsplit; i++) {
203       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
204       ilink = ilink->next;
205     }
206     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
207   }
208   PetscFunctionReturn(0);
209 }
210 
211 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
212     (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
213      VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
214      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
215      VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
216      VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "PCApply_FieldSplit"
220 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
221 {
222   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
223   PetscErrorCode    ierr;
224   PC_FieldSplitLink ilink = jac->head;
225   PetscInt          bs;
226 
227   PetscFunctionBegin;
228   CHKMEMQ;
229   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
230   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
231   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
232 
233   if (jac->type == PC_COMPOSITE_ADDITIVE) {
234     if (jac->defaultsplit) {
235       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
236       while (ilink) {
237 	ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
238 	ilink = ilink->next;
239       }
240       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
241     } else {
242       ierr = VecSet(y,0.0);CHKERRQ(ierr);
243       while (ilink) {
244         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
245 	ilink = ilink->next;
246       }
247     }
248   } else {
249     if (!jac->w1) {
250       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
251       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
252     }
253     ierr = VecSet(y,0.0);CHKERRQ(ierr);
254     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
255     while (ilink->next) {
256       ilink = ilink->next;
257       ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
258       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
259       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
260     }
261     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
262       while (ilink->previous) {
263         ilink = ilink->previous;
264         ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
265         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
266         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
267       }
268     }
269   }
270   CHKMEMQ;
271   PetscFunctionReturn(0);
272 }
273 
274 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
275     (VecScatterBegin(xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
276      VecScatterEnd(xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
277      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
278      VecScatterBegin(ilink->x,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
279      VecScatterEnd(ilink->x,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "PCApply_FieldSplit"
283 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
284 {
285   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
286   PetscErrorCode    ierr;
287   PC_FieldSplitLink ilink = jac->head;
288   PetscInt          bs;
289 
290   PetscFunctionBegin;
291   CHKMEMQ;
292   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
293   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
294   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
295 
296   if (jac->type == PC_COMPOSITE_ADDITIVE) {
297     if (jac->defaultsplit) {
298       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
299       while (ilink) {
300 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
301 	ilink = ilink->next;
302       }
303       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
304     } else {
305       ierr = VecSet(y,0.0);CHKERRQ(ierr);
306       while (ilink) {
307         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
308 	ilink = ilink->next;
309       }
310     }
311   } else {
312     if (!jac->w1) {
313       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
314       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
315     }
316     ierr = VecSet(y,0.0);CHKERRQ(ierr);
317     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
318       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
319       while (ilink->next) {
320         ilink = ilink->next;
321         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
322         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
323         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
324       }
325       while (ilink->previous) {
326         ilink = ilink->previous;
327         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
328         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
329         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
330       }
331     } else {
332       while (ilink->next) {   /* get to last entry in linked list */
333 	ilink = ilink->next;
334       }
335       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
336       while (ilink->previous) {
337 	ilink = ilink->previous;
338 	ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
339 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
340 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
341       }
342     }
343   }
344   CHKMEMQ;
345   PetscFunctionReturn(0);
346 }
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "PCDestroy_FieldSplit"
350 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
351 {
352   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
353   PetscErrorCode    ierr;
354   PC_FieldSplitLink ilink = jac->head,next;
355 
356   PetscFunctionBegin;
357   while (ilink) {
358     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
359     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
360     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
361     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
362     next = ilink->next;
363     ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr);
364     ilink = next;
365   }
366   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
367   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
368   if (jac->is) {
369     PetscInt i;
370     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
371     ierr = PetscFree(jac->is);CHKERRQ(ierr);
372   }
373   if (jac->cis) {
374     PetscInt i;
375     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);}
376     ierr = PetscFree(jac->cis);CHKERRQ(ierr);
377   }
378   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
379   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
380   ierr = PetscFree(jac->csize);CHKERRQ(ierr);
381   ierr = PetscFree(jac);CHKERRQ(ierr);
382   PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
387 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
388 {
389   PetscErrorCode ierr;
390   PetscInt       i = 0,nfields,*fields,bs;
391   PetscTruth     flg;
392   char           optionname[128];
393   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
394 
395   PetscFunctionBegin;
396   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
397   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
398   if (flg) {
399     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
400   }
401   if (jac->bs <= 0) {
402     ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr);
403   }
404   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
405   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
406   while (PETSC_TRUE) {
407     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
408     nfields = jac->bs;
409     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
410     if (!flg) break;
411     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
412     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
413   }
414   ierr = PetscFree(fields);CHKERRQ(ierr);
415   ierr = PetscOptionsTail();CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 
419 /*------------------------------------------------------------------------------------*/
420 
421 EXTERN_C_BEGIN
422 #undef __FUNCT__
423 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
424 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
425 {
426   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
427   PetscErrorCode    ierr;
428   PC_FieldSplitLink ilink,next = jac->head;
429   char              prefix[128];
430   PetscInt          i;
431 
432   PetscFunctionBegin;
433   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
434   for (i=0; i<n; i++) {
435     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
436     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
437   }
438   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr);
439   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
440   ilink->nfields = n;
441   ilink->next    = PETSC_NULL;
442   ierr           = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr);
443   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
444 
445   if (pc->prefix) {
446     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
447   } else {
448     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
449   }
450   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
451 
452   if (!next) {
453     jac->head       = ilink;
454     ilink->previous = PETSC_NULL;
455   } else {
456     while (next->next) {
457       next = next->next;
458     }
459     next->next      = ilink;
460     ilink->previous = next;
461   }
462   jac->nsplits++;
463   PetscFunctionReturn(0);
464 }
465 EXTERN_C_END
466 
467 
468 EXTERN_C_BEGIN
469 #undef __FUNCT__
470 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
471 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
472 {
473   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
474   PetscErrorCode    ierr;
475   PetscInt          cnt = 0;
476   PC_FieldSplitLink ilink = jac->head;
477 
478   PetscFunctionBegin;
479   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
480   while (ilink) {
481     (*subksp)[cnt++] = ilink->ksp;
482     ilink = ilink->next;
483   }
484   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
485   *n = jac->nsplits;
486   PetscFunctionReturn(0);
487 }
488 EXTERN_C_END
489 
490 #undef __FUNCT__
491 #define __FUNCT__ "PCFieldSplitSetFields"
492 /*@
493     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
494 
495     Collective on PC
496 
497     Input Parameters:
498 +   pc  - the preconditioner context
499 .   n - the number of fields in this split
500 .   fields - the fields in this split
501 
502     Level: intermediate
503 
504 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
505 
506 @*/
507 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
508 {
509   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
510 
511   PetscFunctionBegin;
512   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
513   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
514   if (f) {
515     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
516   }
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "PCFieldSplitSetBlockSize"
522 /*@
523     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
524       fieldsplit preconditioner. If not set the matrix block size is used.
525 
526     Collective on PC
527 
528     Input Parameters:
529 +   pc  - the preconditioner context
530 -   bs - the block size
531 
532     Level: intermediate
533 
534 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
535 
536 @*/
537 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
538 {
539   PetscErrorCode ierr,(*f)(PC,PetscInt);
540 
541   PetscFunctionBegin;
542   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
543   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
544   if (f) {
545     ierr = (*f)(pc,bs);CHKERRQ(ierr);
546   }
547   PetscFunctionReturn(0);
548 }
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "PCFieldSplitGetSubKSP"
552 /*@C
553    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
554 
555    Collective on KSP
556 
557    Input Parameter:
558 .  pc - the preconditioner context
559 
560    Output Parameters:
561 +  n - the number of split
562 -  pc - the array of KSP contexts
563 
564    Note:
565    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
566 
567    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
568 
569    Level: advanced
570 
571 .seealso: PCFIELDSPLIT
572 @*/
573 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
574 {
575   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
576 
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
579   PetscValidIntPointer(n,2);
580   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
581   if (f) {
582     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
583   } else {
584     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
585   }
586   PetscFunctionReturn(0);
587 }
588 
589 EXTERN_C_BEGIN
590 #undef __FUNCT__
591 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
592 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
593 {
594   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
595 
596   PetscFunctionBegin;
597   jac->type = type;
598   PetscFunctionReturn(0);
599 }
600 EXTERN_C_END
601 
602 EXTERN_C_BEGIN
603 #undef __FUNCT__
604 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
605 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
606 {
607   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
608 
609   PetscFunctionBegin;
610   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
611   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);
612   jac->bs = bs;
613   PetscFunctionReturn(0);
614 }
615 EXTERN_C_END
616 
617 #undef __FUNCT__
618 #define __FUNCT__ "PCFieldSplitSetType"
619 /*@C
620    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
621 
622    Collective on PC
623 
624    Input Parameter:
625 .  pc - the preconditioner context
626 .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
627 
628    Options Database Key:
629 .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
630 
631    Level: Developer
632 
633 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
634 
635 .seealso: PCCompositeSetType()
636 
637 @*/
638 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
639 {
640   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
641 
642   PetscFunctionBegin;
643   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
644   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
645   if (f) {
646     ierr = (*f)(pc,type);CHKERRQ(ierr);
647   }
648   PetscFunctionReturn(0);
649 }
650 
651 /* -------------------------------------------------------------------------------------*/
652 /*MC
653    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
654                   fields or groups of fields
655 
656 
657      To set options on the solvers for each block append -sub_ to all the PC
658         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
659 
660      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
661          and set the options directly on the resulting KSP object
662 
663    Level: intermediate
664 
665    Options Database Keys:
666 +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
667 .   -pc_splitfield_default - automatically add any fields to additional splits that have not
668                               been supplied explicitly by -pc_splitfield_%d_fields
669 .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
670 -   -pc_splitfield_type <additive,multiplicative>
671 
672    Concepts: physics based preconditioners
673 
674 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
675            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
676 M*/
677 
678 EXTERN_C_BEGIN
679 #undef __FUNCT__
680 #define __FUNCT__ "PCCreate_FieldSplit"
681 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
682 {
683   PetscErrorCode ierr;
684   PC_FieldSplit  *jac;
685 
686   PetscFunctionBegin;
687   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
688   ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr);
689   jac->bs        = -1;
690   jac->nsplits   = 0;
691   jac->type      = PC_COMPOSITE_ADDITIVE;
692 
693   pc->data     = (void*)jac;
694 
695   pc->ops->apply             = PCApply_FieldSplit;
696   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
697   pc->ops->setup             = PCSetUp_FieldSplit;
698   pc->ops->destroy           = PCDestroy_FieldSplit;
699   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
700   pc->ops->view              = PCView_FieldSplit;
701   pc->ops->applyrichardson   = 0;
702 
703   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
704                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
705   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
706                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
707   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
708                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
709   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
710                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
711   PetscFunctionReturn(0);
712 }
713 EXTERN_C_END
714 
715 
716