xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 69a612a9be7fbcc0ab27d96ce41c4e3795ac7f2a)
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 {
8*69a612a9SBarry Smith   KSP               ksp;
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 {
1697bbdb24SBarry Smith   PetscTruth        defaultsplit;
170971522cSBarry Smith   PetscInt          bs;
180971522cSBarry Smith   PetscInt          nsplits;
190971522cSBarry Smith   Vec               *x,*y;
2097bbdb24SBarry Smith   Mat               *pmat;
2197bbdb24SBarry Smith   IS                *is;
2297bbdb24SBarry 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);
39*69a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP 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);
47*69a612a9SBarry Smith       ierr = KSPView(link->ksp,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__
58*69a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
59*69a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(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*69a612a9SBarry Smith   PetscInt          i;
650971522cSBarry Smith 
660971522cSBarry Smith   PetscFunctionBegin;
67*69a612a9SBarry Smith   /* user has not split fields so use default */
68*69a612a9SBarry Smith   if (!link) {
69521d7252SBarry Smith     if (jac->bs <= 0) {
70521d7252SBarry Smith       ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
71521d7252SBarry Smith     }
72521d7252SBarry Smith     for (i=0; i<jac->bs; i++) {
73521d7252SBarry Smith       ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
74521d7252SBarry Smith     }
75521d7252SBarry Smith     link              = jac->head;
7697bbdb24SBarry Smith     jac->defaultsplit = PETSC_TRUE;
77521d7252SBarry Smith   }
78*69a612a9SBarry Smith   PetscFunctionReturn(0);
79*69a612a9SBarry Smith }
80*69a612a9SBarry Smith 
81*69a612a9SBarry Smith 
82*69a612a9SBarry Smith #undef __FUNCT__
83*69a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
84*69a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
85*69a612a9SBarry Smith {
86*69a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
87*69a612a9SBarry Smith   PetscErrorCode    ierr;
88*69a612a9SBarry Smith   PC_FieldSplitLink link;
89*69a612a9SBarry Smith   PetscInt          i,nsplit;
90*69a612a9SBarry Smith   MatStructure      flag = pc->flag;
91*69a612a9SBarry Smith 
92*69a612a9SBarry Smith   PetscFunctionBegin;
93*69a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
9497bbdb24SBarry Smith   nsplit = jac->nsplits;
95*69a612a9SBarry Smith   link   = jac->head;
9697bbdb24SBarry Smith 
9797bbdb24SBarry Smith   /* get the matrices for each split */
9897bbdb24SBarry Smith   if (!jac->is) {
9997bbdb24SBarry Smith     if (jac->defaultsplit) {
10097bbdb24SBarry Smith       PetscInt rstart,rend,bs = nsplit;
10197bbdb24SBarry Smith 
10297bbdb24SBarry Smith       ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
10397bbdb24SBarry Smith       ierr = PetscMalloc(bs*sizeof(IS),&jac->is);CHKERRQ(ierr);
10497bbdb24SBarry Smith       for (i=0; i<bs; i++) {
10597bbdb24SBarry Smith 	ierr = ISCreateStride(pc->comm,(rend-rstart)/bs,rstart+i,bs,&jac->is[i]);CHKERRQ(ierr);
10697bbdb24SBarry Smith       }
10797bbdb24SBarry Smith     } else {
10897bbdb24SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Do not yet support nontrivial split");
10997bbdb24SBarry Smith     }
11097bbdb24SBarry Smith   }
11197bbdb24SBarry Smith   if (!jac->pmat) {
11297bbdb24SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);CHKERRQ(ierr);
11397bbdb24SBarry Smith   } else {
11497bbdb24SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);CHKERRQ(ierr);
11597bbdb24SBarry Smith   }
11697bbdb24SBarry Smith 
11797bbdb24SBarry Smith   /* set up the individual PCs */
11897bbdb24SBarry Smith   i = 0;
1190971522cSBarry Smith   while (link) {
120*69a612a9SBarry Smith     ierr = KSPSetOperators(link->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
121*69a612a9SBarry Smith     ierr = KSPSetFromOptions(link->ksp);CHKERRQ(ierr);
122*69a612a9SBarry Smith     ierr = KSPSetUp(link->ksp);CHKERRQ(ierr);
12397bbdb24SBarry Smith     i++;
1240971522cSBarry Smith     link = link->next;
1250971522cSBarry Smith   }
12697bbdb24SBarry Smith 
12797bbdb24SBarry Smith   /* create work vectors for each split */
12897bbdb24SBarry Smith   if (jac->defaultsplit && !jac->x) {
12997bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
13097bbdb24SBarry Smith     link = jac->head;
13197bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
13297bbdb24SBarry Smith       Mat A;
133*69a612a9SBarry Smith       ierr      = KSPGetOperators(link->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
13497bbdb24SBarry Smith       ierr      = MatGetVecs(A,&link->x,&link->y);CHKERRQ(ierr);
13597bbdb24SBarry Smith       jac->x[i] = link->x;
13697bbdb24SBarry Smith       jac->y[i] = link->y;
13797bbdb24SBarry Smith       link      = link->next;
13897bbdb24SBarry Smith     }
13997bbdb24SBarry Smith   }
1400971522cSBarry Smith   PetscFunctionReturn(0);
1410971522cSBarry Smith }
1420971522cSBarry Smith 
1430971522cSBarry Smith #undef __FUNCT__
1440971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
1450971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
1460971522cSBarry Smith {
1470971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1480971522cSBarry Smith   PetscErrorCode    ierr;
1490971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
1500971522cSBarry Smith 
1510971522cSBarry Smith   PetscFunctionBegin;
1520971522cSBarry Smith   ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
1530971522cSBarry Smith   while (link) {
154*69a612a9SBarry Smith     ierr = KSPSolve(link->ksp,link->x,link->y);CHKERRQ(ierr);
1550971522cSBarry Smith     link = link->next;
1560971522cSBarry Smith   }
1570971522cSBarry Smith   ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
1580971522cSBarry Smith   PetscFunctionReturn(0);
1590971522cSBarry Smith }
1600971522cSBarry Smith 
1610971522cSBarry Smith #undef __FUNCT__
1620971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
1630971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1640971522cSBarry Smith {
1650971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1660971522cSBarry Smith   PetscErrorCode    ierr;
16797bbdb24SBarry Smith   PC_FieldSplitLink link = jac->head,next;
1680971522cSBarry Smith 
1690971522cSBarry Smith   PetscFunctionBegin;
1700971522cSBarry Smith   while (link) {
171*69a612a9SBarry Smith     ierr = KSPDestroy(link->ksp);CHKERRQ(ierr);
172*69a612a9SBarry Smith     if (link->x) {ierr = VecDestroy(link->x);CHKERRQ(ierr);}
173*69a612a9SBarry Smith     if (link->y) {ierr = VecDestroy(link->y);CHKERRQ(ierr);}
17497bbdb24SBarry Smith     next = link->next;
1750971522cSBarry Smith     ierr = PetscFree2(link,link->fields);CHKERRQ(ierr);
17697bbdb24SBarry Smith     link = next;
1770971522cSBarry Smith   }
1780971522cSBarry Smith   if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);}
17997bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
18097bbdb24SBarry Smith   if (jac->is) {
18197bbdb24SBarry Smith     PetscInt i;
18297bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
18397bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
18497bbdb24SBarry Smith   }
1850971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1860971522cSBarry Smith   PetscFunctionReturn(0);
1870971522cSBarry Smith }
1880971522cSBarry Smith 
1890971522cSBarry Smith #undef __FUNCT__
1900971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
1910971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
1920971522cSBarry Smith {
1930971522cSBarry Smith   PetscFunctionBegin;
1940971522cSBarry Smith   PetscFunctionReturn(0);
1950971522cSBarry Smith }
1960971522cSBarry Smith 
1970971522cSBarry Smith /*------------------------------------------------------------------------------------*/
1980971522cSBarry Smith 
1990971522cSBarry Smith EXTERN_C_BEGIN
2000971522cSBarry Smith #undef __FUNCT__
2010971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
2020971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
2030971522cSBarry Smith {
20497bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2050971522cSBarry Smith   PetscErrorCode    ierr;
2060971522cSBarry Smith   PC_FieldSplitLink link,next = jac->head;
207*69a612a9SBarry Smith   char              prefix[128];
2080971522cSBarry Smith 
2090971522cSBarry Smith   PetscFunctionBegin;
2100971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
21197bbdb24SBarry Smith   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&link,n,PetscInt,&link->fields);CHKERRQ(ierr);
21297bbdb24SBarry Smith   ierr = PetscMemcpy(link->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
2130971522cSBarry Smith   link->nfields = n;
21497bbdb24SBarry Smith   link->next    = PETSC_NULL;
215*69a612a9SBarry Smith   ierr          = KSPCreate(pc->comm,&link->ksp);CHKERRQ(ierr);
216*69a612a9SBarry Smith   ierr          = KSPSetType(link->ksp,KSPPREONLY);CHKERRQ(ierr);
217*69a612a9SBarry Smith 
218*69a612a9SBarry Smith   if (pc->prefix) {
219*69a612a9SBarry Smith     sprintf(prefix,"%s_fieldsplit_%d_",pc->prefix,jac->nsplits);
220*69a612a9SBarry Smith   } else {
221*69a612a9SBarry Smith     sprintf(prefix,"fieldsplit_%d_",jac->nsplits);
222*69a612a9SBarry Smith   }
223*69a612a9SBarry Smith   ierr = KSPSetOptionsPrefix(link->ksp,prefix);CHKERRQ(ierr);
2240971522cSBarry Smith 
2250971522cSBarry Smith   if (!next) {
2260971522cSBarry Smith     jac->head = link;
2270971522cSBarry Smith   } else {
2280971522cSBarry Smith     while (next->next) {
2290971522cSBarry Smith       next = next->next;
2300971522cSBarry Smith     }
2310971522cSBarry Smith     next->next = link;
2320971522cSBarry Smith   }
2330971522cSBarry Smith   jac->nsplits++;
2340971522cSBarry Smith   PetscFunctionReturn(0);
2350971522cSBarry Smith }
2360971522cSBarry Smith EXTERN_C_END
2370971522cSBarry Smith 
2380971522cSBarry Smith 
2390971522cSBarry Smith EXTERN_C_BEGIN
2400971522cSBarry Smith #undef __FUNCT__
241*69a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
242*69a612a9SBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
2430971522cSBarry Smith {
2440971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2450971522cSBarry Smith   PetscErrorCode    ierr;
2460971522cSBarry Smith   PetscInt          cnt = 0;
2470971522cSBarry Smith   PC_FieldSplitLink link;
2480971522cSBarry Smith 
2490971522cSBarry Smith   PetscFunctionBegin;
250*69a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
2510971522cSBarry Smith   while (link) {
252*69a612a9SBarry Smith     (*subksp)[cnt++] = link->ksp;
2530971522cSBarry Smith     link = link->next;
2540971522cSBarry Smith   }
2550971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
2560971522cSBarry Smith   *n = jac->nsplits;
2570971522cSBarry Smith   PetscFunctionReturn(0);
2580971522cSBarry Smith }
2590971522cSBarry Smith EXTERN_C_END
2600971522cSBarry Smith 
2610971522cSBarry Smith #undef __FUNCT__
2620971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
2630971522cSBarry Smith /*@
2640971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
2650971522cSBarry Smith 
2660971522cSBarry Smith     Collective on PC
2670971522cSBarry Smith 
2680971522cSBarry Smith     Input Parameters:
2690971522cSBarry Smith +   pc  - the preconditioner context
2700971522cSBarry Smith .   n - the number of fields in this split
2710971522cSBarry Smith .   fields - the fields in this split
2720971522cSBarry Smith 
2730971522cSBarry Smith     Level: intermediate
2740971522cSBarry Smith 
275*69a612a9SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
2760971522cSBarry Smith 
2770971522cSBarry Smith @*/
2780971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
2790971522cSBarry Smith {
2800971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
2810971522cSBarry Smith 
2820971522cSBarry Smith   PetscFunctionBegin;
2830971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
2840971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
2850971522cSBarry Smith   if (f) {
2860971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
2870971522cSBarry Smith   }
2880971522cSBarry Smith   PetscFunctionReturn(0);
2890971522cSBarry Smith }
2900971522cSBarry Smith 
2910971522cSBarry Smith #undef __FUNCT__
292*69a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
2930971522cSBarry Smith /*@C
294*69a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
2950971522cSBarry Smith 
296*69a612a9SBarry Smith    Collective on KSP
2970971522cSBarry Smith 
2980971522cSBarry Smith    Input Parameter:
2990971522cSBarry Smith .  pc - the preconditioner context
3000971522cSBarry Smith 
3010971522cSBarry Smith    Output Parameters:
3020971522cSBarry Smith +  n - the number of split
303*69a612a9SBarry Smith -  pc - the array of KSP contexts
3040971522cSBarry Smith 
3050971522cSBarry Smith    Note:
306*69a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
3070971522cSBarry Smith 
308*69a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
3090971522cSBarry Smith 
3100971522cSBarry Smith    Level: advanced
3110971522cSBarry Smith 
3120971522cSBarry Smith .seealso: PCFIELDSPLIT
3130971522cSBarry Smith @*/
314*69a612a9SBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
3150971522cSBarry Smith {
316*69a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
3170971522cSBarry Smith 
3180971522cSBarry Smith   PetscFunctionBegin;
3190971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3200971522cSBarry Smith   PetscValidIntPointer(n,2);
321*69a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
3220971522cSBarry Smith   if (f) {
323*69a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
3240971522cSBarry Smith   } else {
325*69a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
3260971522cSBarry Smith   }
3270971522cSBarry Smith   PetscFunctionReturn(0);
3280971522cSBarry Smith }
3290971522cSBarry Smith 
3300971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
3310971522cSBarry Smith /*MC
3320971522cSBarry Smith    PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual
3330971522cSBarry Smith                   fields or groups of fields
3340971522cSBarry Smith 
3350971522cSBarry Smith 
3360971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
3370971522cSBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1
3380971522cSBarry Smith 
339*69a612a9SBarry Smith      To set the options on the solvers seperate for each block call PCFieldSplitGetSubKSP()
340*69a612a9SBarry Smith          and set the options directly on the resulting KSP object
3410971522cSBarry Smith 
3420971522cSBarry Smith    Level: intermediate
3430971522cSBarry Smith 
3440971522cSBarry Smith    Concepts: physics based preconditioners
3450971522cSBarry Smith 
3460971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
347*69a612a9SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields()
3480971522cSBarry Smith M*/
3490971522cSBarry Smith 
3500971522cSBarry Smith EXTERN_C_BEGIN
3510971522cSBarry Smith #undef __FUNCT__
3520971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
3530971522cSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc)
3540971522cSBarry Smith {
3550971522cSBarry Smith   PetscErrorCode ierr;
3560971522cSBarry Smith   PC_FieldSplit  *jac;
3570971522cSBarry Smith 
3580971522cSBarry Smith   PetscFunctionBegin;
3590971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
3600971522cSBarry Smith   PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));
3610971522cSBarry Smith   ierr = PetscMemzero(jac,sizeof(PC_FieldSplit));CHKERRQ(ierr);
3620971522cSBarry Smith   jac->bs      = -1;
3630971522cSBarry Smith   jac->nsplits = 0;
3640971522cSBarry Smith   pc->data     = (void*)jac;
3650971522cSBarry Smith 
3660971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
3670971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
3680971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
3690971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
3700971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
3710971522cSBarry Smith   pc->ops->applyrichardson   = 0;
3720971522cSBarry Smith 
373*69a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
374*69a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
3750971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
3760971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
3770971522cSBarry Smith   PetscFunctionReturn(0);
3780971522cSBarry Smith }
3790971522cSBarry Smith EXTERN_C_END
3800971522cSBarry Smith 
3810971522cSBarry Smith 
382