xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 97bbdb249ec3ac0a9d48a4883afe096645382371)
10971522cSBarry Smith /*
20971522cSBarry Smith 
30971522cSBarry Smith */
40971522cSBarry Smith #include "src/ksp/pc/pcimpl.h"     /*I "petscpc.h" I*/
50971522cSBarry Smith 
60971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
70971522cSBarry Smith struct _PC_FieldSplitLink {
80971522cSBarry Smith   PC                pc;
90971522cSBarry Smith   Vec               x,y;
100971522cSBarry Smith   PetscInt          nfields;
110971522cSBarry Smith   PetscInt          *fields;
120971522cSBarry Smith   PC_FieldSplitLink next;
130971522cSBarry Smith };
140971522cSBarry Smith 
150971522cSBarry Smith typedef struct {
16*97bbdb24SBarry Smith   PetscTruth        defaultsplit;
170971522cSBarry Smith   PetscInt          bs;
180971522cSBarry Smith   PetscInt          nsplits;
190971522cSBarry Smith   Vec               *x,*y;
20*97bbdb24SBarry Smith   Mat               *pmat;
21*97bbdb24SBarry Smith   IS                *is;
22*97bbdb24SBarry Smith   PC_FieldSplitLink head;
230971522cSBarry Smith } PC_FieldSplit;
240971522cSBarry Smith 
250971522cSBarry Smith #undef __FUNCT__
260971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
270971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
280971522cSBarry Smith {
290971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
300971522cSBarry Smith   PetscErrorCode    ierr;
310971522cSBarry Smith   PetscTruth        iascii;
320971522cSBarry Smith   PetscInt          i,j;
330971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
340971522cSBarry Smith 
350971522cSBarry Smith   PetscFunctionBegin;
360971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
370971522cSBarry Smith   if (iascii) {
380971522cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit: total splits = %D",jac->nsplits);CHKERRQ(ierr);
390971522cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following PC objects:\n");CHKERRQ(ierr);
400971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
410971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
420971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
430971522cSBarry Smith       for (j=0; j<link->nfields; j++) {
440971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%D \n",link->fields[j]);CHKERRQ(ierr);
450971522cSBarry Smith       }
460971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
470971522cSBarry Smith       ierr = PCView(link->pc,viewer);CHKERRQ(ierr);
480971522cSBarry Smith       link = link->next;
490971522cSBarry Smith     }
500971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
510971522cSBarry Smith   } else {
520971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
530971522cSBarry Smith   }
540971522cSBarry Smith   PetscFunctionReturn(0);
550971522cSBarry Smith }
560971522cSBarry Smith 
570971522cSBarry Smith #undef __FUNCT__
580971522cSBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
590971522cSBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
600971522cSBarry Smith {
610971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
620971522cSBarry Smith   PetscErrorCode    ierr;
630971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
64*97bbdb24SBarry Smith   PetscInt          i,nsplit;
65*97bbdb24SBarry Smith   MatStructure      flag = pc->flag;
660971522cSBarry Smith 
670971522cSBarry Smith   PetscFunctionBegin;
68521d7252SBarry Smith   if (jac->bs <= 0) {
69521d7252SBarry Smith     ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
70521d7252SBarry Smith   }
71*97bbdb24SBarry Smith 
72*97bbdb24SBarry Smith   /* user has not split fields so use default */
73*97bbdb24SBarry Smith   if (!link) {
74521d7252SBarry Smith     for (i=0; i<jac->bs; i++) {
75521d7252SBarry Smith       ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
76521d7252SBarry Smith     }
77521d7252SBarry Smith     link              = jac->head;
78*97bbdb24SBarry Smith     jac->defaultsplit = PETSC_TRUE;
79521d7252SBarry Smith   }
80*97bbdb24SBarry Smith   nsplit = jac->nsplits;
81*97bbdb24SBarry Smith 
82*97bbdb24SBarry Smith   /* get the matrices for each split */
83*97bbdb24SBarry Smith   if (!jac->is) {
84*97bbdb24SBarry Smith     if (jac->defaultsplit) {
85*97bbdb24SBarry Smith       PetscInt rstart,rend,bs = nsplit;
86*97bbdb24SBarry Smith 
87*97bbdb24SBarry Smith       ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
88*97bbdb24SBarry Smith       ierr = PetscMalloc(bs*sizeof(IS),&jac->is);CHKERRQ(ierr);
89*97bbdb24SBarry Smith       for (i=0; i<bs; i++) {
90*97bbdb24SBarry Smith 	ierr = ISCreateStride(pc->comm,(rend-rstart)/bs,rstart+i,bs,&jac->is[i]);CHKERRQ(ierr);
91*97bbdb24SBarry Smith       }
92*97bbdb24SBarry Smith     } else {
93*97bbdb24SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Do not yet support nontrivial split");
94*97bbdb24SBarry Smith     }
95*97bbdb24SBarry Smith   }
96*97bbdb24SBarry Smith   if (!jac->pmat) {
97*97bbdb24SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);CHKERRQ(ierr);
98*97bbdb24SBarry Smith   } else {
99*97bbdb24SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);CHKERRQ(ierr);
100*97bbdb24SBarry Smith   }
101*97bbdb24SBarry Smith 
102*97bbdb24SBarry Smith   /* set up the individual PCs */
103*97bbdb24SBarry Smith   i = 0;
1040971522cSBarry Smith   while (link) {
105*97bbdb24SBarry Smith     ierr = PCSetOperators(link->pc,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
1060971522cSBarry Smith     ierr = PCSetUp(link->pc);CHKERRQ(ierr);
107*97bbdb24SBarry Smith     i++;
1080971522cSBarry Smith     link = link->next;
1090971522cSBarry Smith   }
110*97bbdb24SBarry Smith 
111*97bbdb24SBarry Smith   /* create work vectors for each split */
112*97bbdb24SBarry Smith   if (jac->defaultsplit && !jac->x) {
113*97bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
114*97bbdb24SBarry Smith     link = jac->head;
115*97bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
116*97bbdb24SBarry Smith       Mat A;
117*97bbdb24SBarry Smith       ierr      = PCGetOperators(link->pc,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
118*97bbdb24SBarry Smith       ierr      = MatGetVecs(A,&link->x,&link->y);CHKERRQ(ierr);
119*97bbdb24SBarry Smith       jac->x[i] = link->x;
120*97bbdb24SBarry Smith       jac->y[i] = link->y;
121*97bbdb24SBarry Smith       link      = link->next;
122*97bbdb24SBarry Smith     }
123*97bbdb24SBarry Smith   }
1240971522cSBarry Smith   PetscFunctionReturn(0);
1250971522cSBarry Smith }
1260971522cSBarry Smith 
1270971522cSBarry Smith #undef __FUNCT__
1280971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
1290971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
1300971522cSBarry Smith {
1310971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1320971522cSBarry Smith   PetscErrorCode    ierr;
1330971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
1340971522cSBarry Smith 
1350971522cSBarry Smith   PetscFunctionBegin;
1360971522cSBarry Smith   ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
1370971522cSBarry Smith   while (link) {
1380971522cSBarry Smith     ierr = PCApply(link->pc,link->x,link->y);CHKERRQ(ierr);
1390971522cSBarry Smith     link = link->next;
1400971522cSBarry Smith   }
1410971522cSBarry Smith   ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
1420971522cSBarry Smith   PetscFunctionReturn(0);
1430971522cSBarry Smith }
1440971522cSBarry Smith 
1450971522cSBarry Smith #undef __FUNCT__
1460971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
1470971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1480971522cSBarry Smith {
1490971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1500971522cSBarry Smith   PetscErrorCode    ierr;
151*97bbdb24SBarry Smith   PC_FieldSplitLink link = jac->head,next;
1520971522cSBarry Smith 
1530971522cSBarry Smith   PetscFunctionBegin;
1540971522cSBarry Smith   while (link) {
1550971522cSBarry Smith     ierr = PCDestroy(link->pc);CHKERRQ(ierr);
156*97bbdb24SBarry Smith     ierr = VecDestroy(link->x);CHKERRQ(ierr);
157*97bbdb24SBarry Smith     ierr = VecDestroy(link->y);CHKERRQ(ierr);
158*97bbdb24SBarry Smith     next = link->next;
1590971522cSBarry Smith     ierr = PetscFree2(link,link->fields);CHKERRQ(ierr);
160*97bbdb24SBarry Smith     link = next;
1610971522cSBarry Smith   }
1620971522cSBarry Smith   if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);}
163*97bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
164*97bbdb24SBarry Smith   if (jac->is) {
165*97bbdb24SBarry Smith     PetscInt i;
166*97bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
167*97bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
168*97bbdb24SBarry Smith   }
1690971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1700971522cSBarry Smith   PetscFunctionReturn(0);
1710971522cSBarry Smith }
1720971522cSBarry Smith 
1730971522cSBarry Smith #undef __FUNCT__
1740971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
1750971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
1760971522cSBarry Smith {
1770971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1780971522cSBarry Smith   PetscErrorCode    ierr;
1790971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
1800971522cSBarry Smith 
1810971522cSBarry Smith   PetscFunctionBegin;
1820971522cSBarry Smith   while (link) {
183*97bbdb24SBarry Smith     ierr = PCSetFromOptions(link->pc);CHKERRQ(ierr);
1840971522cSBarry Smith     link = link->next;
1850971522cSBarry Smith   }
1860971522cSBarry Smith   PetscFunctionReturn(0);
1870971522cSBarry Smith }
1880971522cSBarry Smith 
1890971522cSBarry Smith /*------------------------------------------------------------------------------------*/
1900971522cSBarry Smith 
1910971522cSBarry Smith EXTERN_C_BEGIN
1920971522cSBarry Smith #undef __FUNCT__
1930971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
1940971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
1950971522cSBarry Smith {
196*97bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1970971522cSBarry Smith   PetscErrorCode    ierr;
1980971522cSBarry Smith   PC_FieldSplitLink link,next = jac->head;
1990971522cSBarry Smith 
2000971522cSBarry Smith   PetscFunctionBegin;
2010971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
202*97bbdb24SBarry Smith   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&link,n,PetscInt,&link->fields);CHKERRQ(ierr);
203*97bbdb24SBarry Smith   ierr = PetscMemcpy(link->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
2040971522cSBarry Smith   link->nfields = n;
205*97bbdb24SBarry Smith   link->next    = PETSC_NULL;
206*97bbdb24SBarry Smith   ierr          = PCCreate(pc->comm,&link->pc);CHKERRQ(ierr);
2070971522cSBarry Smith 
2080971522cSBarry Smith   if (!next) {
2090971522cSBarry Smith     jac->head = link;
2100971522cSBarry Smith   } else {
2110971522cSBarry Smith     while (next->next) {
2120971522cSBarry Smith       next = next->next;
2130971522cSBarry Smith     }
2140971522cSBarry Smith     next->next = link;
2150971522cSBarry Smith   }
2160971522cSBarry Smith   jac->nsplits++;
2170971522cSBarry Smith   PetscFunctionReturn(0);
2180971522cSBarry Smith }
2190971522cSBarry Smith EXTERN_C_END
2200971522cSBarry Smith 
2210971522cSBarry Smith 
2220971522cSBarry Smith EXTERN_C_BEGIN
2230971522cSBarry Smith #undef __FUNCT__
2240971522cSBarry Smith #define __FUNCT__ "PCFieldSplitGetSubPC_FieldSplit"
2250971522cSBarry Smith PetscErrorCode PCFieldSplitGetSubPC_FieldSplit(PC pc,PetscInt *n,PC **subpc)
2260971522cSBarry Smith {
2270971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2280971522cSBarry Smith   PetscErrorCode    ierr;
2290971522cSBarry Smith   PetscInt          cnt = 0;
2300971522cSBarry Smith   PC_FieldSplitLink link;
2310971522cSBarry Smith 
2320971522cSBarry Smith   PetscFunctionBegin;
2330971522cSBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(PC*),subpc);CHKERRQ(ierr);
2340971522cSBarry Smith   while (link) {
2350971522cSBarry Smith     (*subpc)[cnt++] = link->pc;
2360971522cSBarry Smith     link = link->next;
2370971522cSBarry Smith   }
2380971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
2390971522cSBarry Smith   *n = jac->nsplits;
2400971522cSBarry Smith   PetscFunctionReturn(0);
2410971522cSBarry Smith }
2420971522cSBarry Smith EXTERN_C_END
2430971522cSBarry Smith 
2440971522cSBarry Smith #undef __FUNCT__
2450971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
2460971522cSBarry Smith /*@
2470971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
2480971522cSBarry Smith 
2490971522cSBarry Smith     Collective on PC
2500971522cSBarry Smith 
2510971522cSBarry Smith     Input Parameters:
2520971522cSBarry Smith +   pc  - the preconditioner context
2530971522cSBarry Smith .   n - the number of fields in this split
2540971522cSBarry Smith .   fields - the fields in this split
2550971522cSBarry Smith 
2560971522cSBarry Smith     Level: intermediate
2570971522cSBarry Smith 
2583d30b1ffSBarry Smith .seealso: PCFieldSplitGetSubPC(), PCFIELDSPLIT
2590971522cSBarry Smith 
2600971522cSBarry Smith @*/
2610971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
2620971522cSBarry Smith {
2630971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
2640971522cSBarry Smith 
2650971522cSBarry Smith   PetscFunctionBegin;
2660971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
2670971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
2680971522cSBarry Smith   if (f) {
2690971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
2700971522cSBarry Smith   }
2710971522cSBarry Smith   PetscFunctionReturn(0);
2720971522cSBarry Smith }
2730971522cSBarry Smith 
2740971522cSBarry Smith #undef __FUNCT__
2750971522cSBarry Smith #define __FUNCT__ "PCFieldSplitGetSubPC"
2760971522cSBarry Smith /*@C
2770971522cSBarry Smith    PCFieldSplitGetSubPC - Gets the PC contexts for all splits
2780971522cSBarry Smith 
2790971522cSBarry Smith    Collective on PC
2800971522cSBarry Smith 
2810971522cSBarry Smith    Input Parameter:
2820971522cSBarry Smith .  pc - the preconditioner context
2830971522cSBarry Smith 
2840971522cSBarry Smith    Output Parameters:
2850971522cSBarry Smith +  n - the number of split
2860971522cSBarry Smith -  pc - the array of PC contexts
2870971522cSBarry Smith 
2880971522cSBarry Smith    Note:
2890971522cSBarry Smith    After PCFieldSplitGetSubPC() the array of PCs IS to be freed
2900971522cSBarry Smith 
2910971522cSBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubPC().
2920971522cSBarry Smith 
2930971522cSBarry Smith    Level: advanced
2940971522cSBarry Smith 
2950971522cSBarry Smith .keywords: PC,
2960971522cSBarry Smith 
2970971522cSBarry Smith .seealso: PCFIELDSPLIT
2980971522cSBarry Smith @*/
2990971522cSBarry Smith PetscErrorCode PCFieldSplitGetSubPC(PC pc,PetscInt *n,PC *subpc[])
3000971522cSBarry Smith {
3010971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,PC **);
3020971522cSBarry Smith 
3030971522cSBarry Smith   PetscFunctionBegin;
3040971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3050971522cSBarry Smith   PetscValidIntPointer(n,2);
3060971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubPC_C",(void (**)(void))&f);CHKERRQ(ierr);
3070971522cSBarry Smith   if (f) {
3080971522cSBarry Smith     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
3090971522cSBarry Smith   } else {
3100971522cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subpc for this type of PC");
3110971522cSBarry Smith   }
3120971522cSBarry Smith   PetscFunctionReturn(0);
3130971522cSBarry Smith }
3140971522cSBarry Smith 
3150971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
3160971522cSBarry Smith /*MC
3170971522cSBarry Smith    PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual
3180971522cSBarry Smith                   fields or groups of fields
3190971522cSBarry Smith 
3200971522cSBarry Smith 
3210971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
3220971522cSBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1
3230971522cSBarry Smith 
3240971522cSBarry Smith      To set the options on the solvers seperate for each block call PCFieldSplitGetSubPC()
3250971522cSBarry Smith          and set the options directly on the resulting PC object
3260971522cSBarry Smith 
3270971522cSBarry Smith    Level: intermediate
3280971522cSBarry Smith 
3290971522cSBarry Smith    Concepts: physics based preconditioners
3300971522cSBarry Smith 
3310971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
3320971522cSBarry Smith            PCFieldSplitGetSubPC(), PCFieldSplitSetFields()
3330971522cSBarry Smith M*/
3340971522cSBarry Smith 
3350971522cSBarry Smith EXTERN_C_BEGIN
3360971522cSBarry Smith #undef __FUNCT__
3370971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
3380971522cSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc)
3390971522cSBarry Smith {
3400971522cSBarry Smith   PetscErrorCode ierr;
3410971522cSBarry Smith   PC_FieldSplit  *jac;
3420971522cSBarry Smith 
3430971522cSBarry Smith   PetscFunctionBegin;
3440971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
3450971522cSBarry Smith   PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));
3460971522cSBarry Smith   ierr = PetscMemzero(jac,sizeof(PC_FieldSplit));CHKERRQ(ierr);
3470971522cSBarry Smith   jac->bs      = -1;
3480971522cSBarry Smith   jac->nsplits = 0;
3490971522cSBarry Smith   pc->data     = (void*)jac;
3500971522cSBarry Smith 
3510971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
3520971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
3530971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
3540971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
3550971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
3560971522cSBarry Smith   pc->ops->applyrichardson   = 0;
3570971522cSBarry Smith 
3580971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubPC_C","PCFieldSplitGetSubPC_FieldSplit",
3590971522cSBarry Smith                     PCFieldSplitGetSubPC_FieldSplit);CHKERRQ(ierr);
3600971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
3610971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
3620971522cSBarry Smith   PetscFunctionReturn(0);
3630971522cSBarry Smith }
3640971522cSBarry Smith EXTERN_C_END
3650971522cSBarry Smith 
3660971522cSBarry Smith 
367