xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 2a6744eb01855f5aa328eb8fdf4b0d01e72ad151)
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;
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\n",PCCompositeTypes[jac->type],jac->nsplits);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 
114   PetscFunctionBegin;
115   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
116   nsplit = jac->nsplits;
117   ilink  = jac->head;
118 
119   /* get the matrices for each split */
120   if (!jac->is) {
121     PetscInt rstart,rend,nslots,bs;
122 
123     ierr   = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
124     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
125     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
126     nslots = (rend - rstart)/bs;
127     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr);
128     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->cis);CHKERRQ(ierr);
129     ierr   = PetscMalloc(nsplit*sizeof(PetscInt),&jac->csize);CHKERRQ(ierr);
130     for (i=0; i<nsplit; i++) {
131       if (jac->defaultsplit) {
132 	ierr     = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr);
133         if (bs != nsplit) SETERRQ2(PETSC_ERR_PLIB,"With default-split the number of fields %D must equal the matrix block size %D",nsplit,bs);
134         jac->csize[i] = ccsize/nsplit;
135       } else {
136         PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
137         PetscTruth sorted;
138         ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
139         for (j=0; j<nslots; j++) {
140           for (k=0; k<nfields; k++) {
141             ii[nfields*j + k] = rstart + bs*j + fields[k];
142           }
143         }
144 	ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr);
145         jac->csize[i] = (ccsize/bs)*ilink->nfields;
146         ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr);
147         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
148         ierr = PetscFree(ii);CHKERRQ(ierr);
149         ilink = ilink->next;
150       }
151       ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr);
152     }
153   }
154 
155   if (!jac->pmat) {
156     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
157     for (i=0; i<nsplit; i++) {
158       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
159     }
160   } else {
161     for (i=0; i<nsplit; i++) {
162       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
163     }
164   }
165 
166   /* set up the individual PCs */
167   i    = 0;
168   ilink = jac->head;
169   while (ilink) {
170     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
171     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
172     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
173     i++;
174     ilink = ilink->next;
175   }
176 
177   /* create work vectors for each split */
178   if (!jac->x) {
179     Vec xtmp;
180     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
181     ilink = jac->head;
182     for (i=0; i<nsplit; i++) {
183       Mat A;
184       ierr      = KSPGetOperators(ilink->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
185       ierr      = MatGetVecs(A,&ilink->x,&ilink->y);CHKERRQ(ierr);
186       jac->x[i] = ilink->x;
187       jac->y[i] = ilink->y;
188       ilink      = ilink->next;
189     }
190     /* compute scatter contexts needed by multiplicative versions and non-default splits */
191 
192     ilink = jac->head;
193     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
194     for (i=0; i<nsplit; i++) {
195       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
196       ilink = ilink->next;
197     }
198     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
199   }
200   PetscFunctionReturn(0);
201 }
202 
203 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
204     (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
205      VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
206      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
207      VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
208      VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "PCApply_FieldSplit"
212 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
213 {
214   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
215   PetscErrorCode    ierr;
216   PC_FieldSplitLink ilink = jac->head;
217   PetscInt          bs;
218 
219   PetscFunctionBegin;
220   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
221   if (bs != jac->bs) {
222     SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector blocksize %D does not match Fieldsplit blocksize %D",bs,jac->bs);
223   }
224   if (jac->type == PC_COMPOSITE_ADDITIVE) {
225     if (jac->defaultsplit) {
226       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
227       while (ilink) {
228 	ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
229 	ilink = ilink->next;
230       }
231       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
232     } else {
233       PetscInt    i = 0;
234 
235       ierr = VecSet(y,0.0);CHKERRQ(ierr);
236       while (ilink) {
237         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
238 	ilink = ilink->next;
239 	i++;
240       }
241     }
242   } else {
243     if (!jac->w1) {
244       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
245       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
246     }
247     ierr = VecSet(y,0.0);CHKERRQ(ierr);
248     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
249     while (ilink->next) {
250       ilink = ilink->next;
251       ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
252       ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
253       ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
254     }
255   }
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "PCDestroy_FieldSplit"
261 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
262 {
263   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
264   PetscErrorCode    ierr;
265   PC_FieldSplitLink ilink = jac->head,next;
266 
267   PetscFunctionBegin;
268   while (ilink) {
269     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
270     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
271     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
272     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
273     next = ilink->next;
274     ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr);
275     ilink = next;
276   }
277   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
278   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
279   if (jac->is) {
280     PetscInt i;
281     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
282     ierr = PetscFree(jac->is);CHKERRQ(ierr);
283   }
284   if (jac->cis) {
285     PetscInt i;
286     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);}
287     ierr = PetscFree(jac->cis);CHKERRQ(ierr);
288   }
289   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
290   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
291   ierr = PetscFree(jac->csize);CHKERRQ(ierr);
292   ierr = PetscFree(jac);CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 
296 #undef __FUNCT__
297 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
298 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
299 {
300   PetscErrorCode ierr;
301   PetscInt       i = 0,nfields,fields[12];
302   PetscTruth     flg;
303   char           optionname[128];
304   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
305 
306   PetscFunctionBegin;
307   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
308   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
309   while (PETSC_TRUE) {
310     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
311     nfields = 12;
312     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
313     if (!flg) break;
314     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
315     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
316   }
317   ierr = PetscOptionsTail();CHKERRQ(ierr);
318   PetscFunctionReturn(0);
319 }
320 
321 /*------------------------------------------------------------------------------------*/
322 
323 EXTERN_C_BEGIN
324 #undef __FUNCT__
325 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
326 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
327 {
328   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
329   PetscErrorCode    ierr;
330   PC_FieldSplitLink ilink,next = jac->head;
331   char              prefix[128];
332 
333   PetscFunctionBegin;
334   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
335   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr);
336   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
337   ilink->nfields = n;
338   ilink->next    = PETSC_NULL;
339   ierr          = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr);
340   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
341 
342   if (pc->prefix) {
343     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
344   } else {
345     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
346   }
347   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
348 
349   if (!next) {
350     jac->head = ilink;
351   } else {
352     while (next->next) {
353       next = next->next;
354     }
355     next->next = ilink;
356   }
357   jac->nsplits++;
358   PetscFunctionReturn(0);
359 }
360 EXTERN_C_END
361 
362 
363 EXTERN_C_BEGIN
364 #undef __FUNCT__
365 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
366 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
367 {
368   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
369   PetscErrorCode    ierr;
370   PetscInt          cnt = 0;
371   PC_FieldSplitLink ilink = jac->head;
372 
373   PetscFunctionBegin;
374   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
375   while (ilink) {
376     (*subksp)[cnt++] = ilink->ksp;
377     ilink = ilink->next;
378   }
379   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
380   *n = jac->nsplits;
381   PetscFunctionReturn(0);
382 }
383 EXTERN_C_END
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "PCFieldSplitSetFields"
387 /*@
388     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
389 
390     Collective on PC
391 
392     Input Parameters:
393 +   pc  - the preconditioner context
394 .   n - the number of fields in this split
395 .   fields - the fields in this split
396 
397     Level: intermediate
398 
399 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
400 
401 @*/
402 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
403 {
404   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
405 
406   PetscFunctionBegin;
407   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
408   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
409   if (f) {
410     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
411   }
412   PetscFunctionReturn(0);
413 }
414 
415 #undef __FUNCT__
416 #define __FUNCT__ "PCFieldSplitGetSubKSP"
417 /*@C
418    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
419 
420    Collective on KSP
421 
422    Input Parameter:
423 .  pc - the preconditioner context
424 
425    Output Parameters:
426 +  n - the number of split
427 -  pc - the array of KSP contexts
428 
429    Note:
430    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
431 
432    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
433 
434    Level: advanced
435 
436 .seealso: PCFIELDSPLIT
437 @*/
438 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
439 {
440   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
441 
442   PetscFunctionBegin;
443   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
444   PetscValidIntPointer(n,2);
445   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
446   if (f) {
447     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
448   } else {
449     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
450   }
451   PetscFunctionReturn(0);
452 }
453 
454 EXTERN_C_BEGIN
455 #undef __FUNCT__
456 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
457 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
458 {
459   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
460 
461   PetscFunctionBegin;
462   jac->type = type;
463   PetscFunctionReturn(0);
464 }
465 EXTERN_C_END
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "PCFieldSplitSetType"
469 /*@C
470    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
471 
472    Collective on PC
473 
474    Input Parameter:
475 .  pc - the preconditioner context
476 .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
477 
478    Options Database Key:
479 .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
480 
481    Level: Developer
482 
483 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
484 
485 .seealso: PCCompositeSetType()
486 
487 @*/
488 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
489 {
490   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
491 
492   PetscFunctionBegin;
493   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
494   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
495   if (f) {
496     ierr = (*f)(pc,type);CHKERRQ(ierr);
497   }
498   PetscFunctionReturn(0);
499 }
500 
501 /* -------------------------------------------------------------------------------------*/
502 /*MC
503    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
504                   fields or groups of fields
505 
506 
507      To set options on the solvers for each block append -sub_ to all the PC
508         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
509 
510      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
511          and set the options directly on the resulting KSP object
512 
513    Level: intermediate
514 
515    Options Database Keys:
516 +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
517 .   -pc_splitfield_default - automatically add any fields to additional splits that have not
518                               been supplied explicitly by -pc_splitfield_%d_fields
519 -   -pc_splitfield_type <additive,multiplicative>
520 
521    Concepts: physics based preconditioners
522 
523 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
524            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
525 M*/
526 
527 EXTERN_C_BEGIN
528 #undef __FUNCT__
529 #define __FUNCT__ "PCCreate_FieldSplit"
530 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
531 {
532   PetscErrorCode ierr;
533   PC_FieldSplit  *jac;
534 
535   PetscFunctionBegin;
536   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
537   ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr);
538   jac->bs      = -1;
539   jac->nsplits = 0;
540   jac->type    = PC_COMPOSITE_ADDITIVE;
541   pc->data     = (void*)jac;
542 
543   pc->ops->apply             = PCApply_FieldSplit;
544   pc->ops->setup             = PCSetUp_FieldSplit;
545   pc->ops->destroy           = PCDestroy_FieldSplit;
546   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
547   pc->ops->view              = PCView_FieldSplit;
548   pc->ops->applyrichardson   = 0;
549 
550   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
551                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
552   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
553                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
554   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
555                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
556   PetscFunctionReturn(0);
557 }
558 EXTERN_C_END
559 
560 
561