xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision cf50294291ee912251d43f98e0337904404b8dbf)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
30971522cSBarry Smith /*
40971522cSBarry Smith 
50971522cSBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
70971522cSBarry Smith 
80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
90971522cSBarry Smith struct _PC_FieldSplitLink {
1069a612a9SBarry Smith   KSP               ksp;
110971522cSBarry Smith   Vec               x,y;
120971522cSBarry Smith   PetscInt          nfields;
130971522cSBarry Smith   PetscInt          *fields;
141b9fc7fcSBarry Smith   VecScatter        sctx;
150971522cSBarry Smith   PC_FieldSplitLink next;
160971522cSBarry Smith };
170971522cSBarry Smith 
180971522cSBarry Smith typedef struct {
1979416396SBarry Smith   PCCompositeType   type;              /* additive or multiplicative */
2097bbdb24SBarry Smith   PetscTruth        defaultsplit;
210971522cSBarry Smith   PetscInt          bs;
22*cf502942SBarry Smith   PetscInt          nsplits,*csize;
2379416396SBarry Smith   Vec               *x,*y,w1,w2;
2497bbdb24SBarry Smith   Mat               *pmat;
25*cf502942SBarry Smith   IS                *is,*cis;
2697bbdb24SBarry Smith   PC_FieldSplitLink head;
270971522cSBarry Smith } PC_FieldSplit;
280971522cSBarry Smith 
290971522cSBarry Smith #undef __FUNCT__
300971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
310971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
320971522cSBarry Smith {
330971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
340971522cSBarry Smith   PetscErrorCode    ierr;
350971522cSBarry Smith   PetscTruth        iascii;
360971522cSBarry Smith   PetscInt          i,j;
375a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
380971522cSBarry Smith 
390971522cSBarry Smith   PetscFunctionBegin;
400971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
410971522cSBarry Smith   if (iascii) {
42*cf502942SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
4369a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
440971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
450971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
460971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
4779416396SBarry Smith       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
485a9f2f41SSatish Balay       for (j=0; j<ilink->nfields; j++) {
4979416396SBarry Smith         if (j > 0) {
5079416396SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
5179416396SBarry Smith         }
525a9f2f41SSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
530971522cSBarry Smith       }
540971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
5579416396SBarry Smith       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
565a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
575a9f2f41SSatish Balay       ilink = ilink->next;
580971522cSBarry Smith     }
590971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
600971522cSBarry Smith   } else {
610971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
620971522cSBarry Smith   }
630971522cSBarry Smith   PetscFunctionReturn(0);
640971522cSBarry Smith }
650971522cSBarry Smith 
660971522cSBarry Smith #undef __FUNCT__
6769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
6869a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
690971522cSBarry Smith {
700971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
710971522cSBarry Smith   PetscErrorCode    ierr;
725a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
7369a612a9SBarry Smith   PetscInt          i;
7479416396SBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields;
750971522cSBarry Smith 
760971522cSBarry Smith   PetscFunctionBegin;
774dc4c822SBarry Smith   ierr = PetscOptionsGetTruth(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
785a9f2f41SSatish Balay   if (!ilink || flg) {
79ae15b995SBarry Smith     ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
80521d7252SBarry Smith     if (jac->bs <= 0) {
81521d7252SBarry Smith       ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
82521d7252SBarry Smith     }
8379416396SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
8479416396SBarry Smith     ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
855a9f2f41SSatish Balay     while (ilink) {
865a9f2f41SSatish Balay       for (i=0; i<ilink->nfields; i++) {
875a9f2f41SSatish Balay         fields[ilink->fields[i]] = PETSC_TRUE;
88521d7252SBarry Smith       }
895a9f2f41SSatish Balay       ilink = ilink->next;
9079416396SBarry Smith     }
9197bbdb24SBarry Smith     jac->defaultsplit = PETSC_TRUE;
9279416396SBarry Smith     for (i=0; i<jac->bs; i++) {
9379416396SBarry Smith       if (!fields[i]) {
9479416396SBarry Smith 	ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
9579416396SBarry Smith       } else {
9679416396SBarry Smith         jac->defaultsplit = PETSC_FALSE;
9779416396SBarry Smith       }
9879416396SBarry Smith     }
99521d7252SBarry Smith   }
10069a612a9SBarry Smith   PetscFunctionReturn(0);
10169a612a9SBarry Smith }
10269a612a9SBarry Smith 
10369a612a9SBarry Smith 
10469a612a9SBarry Smith #undef __FUNCT__
10569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
10669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
10769a612a9SBarry Smith {
10869a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10969a612a9SBarry Smith   PetscErrorCode    ierr;
1105a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
111*cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
11269a612a9SBarry Smith   MatStructure      flag = pc->flag;
11369a612a9SBarry Smith 
11469a612a9SBarry Smith   PetscFunctionBegin;
11569a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
11697bbdb24SBarry Smith   nsplit = jac->nsplits;
1175a9f2f41SSatish Balay   ilink  = jac->head;
11897bbdb24SBarry Smith 
11997bbdb24SBarry Smith   /* get the matrices for each split */
12097bbdb24SBarry Smith   if (!jac->is) {
1211b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
12297bbdb24SBarry Smith 
1231b9fc7fcSBarry Smith     ierr   = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
12497bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
125*cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
1261b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1271b9fc7fcSBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr);
128*cf502942SBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->cis);CHKERRQ(ierr);
129*cf502942SBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(PetscInt),&jac->csize);CHKERRQ(ierr);
1301b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1311b9fc7fcSBarry Smith       if (jac->defaultsplit) {
1321b9fc7fcSBarry Smith 	ierr     = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr);
133*cf502942SBarry Smith         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*cf502942SBarry Smith         jac->csize[i] = ccsize/nsplit;
13597bbdb24SBarry Smith       } else {
1365a9f2f41SSatish Balay         PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1371b9fc7fcSBarry Smith         PetscTruth sorted;
1385a9f2f41SSatish Balay         ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1391b9fc7fcSBarry Smith         for (j=0; j<nslots; j++) {
1401b9fc7fcSBarry Smith           for (k=0; k<nfields; k++) {
1411b9fc7fcSBarry Smith             ii[nfields*j + k] = rstart + bs*j + fields[k];
14297bbdb24SBarry Smith           }
14397bbdb24SBarry Smith         }
1441b9fc7fcSBarry Smith 	ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr);
145*cf502942SBarry Smith         jac->csize[i] = (ccsize/bs)*ilink->nfields;
1461b9fc7fcSBarry Smith         ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr);
1471b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1481b9fc7fcSBarry Smith         ierr = PetscFree(ii);CHKERRQ(ierr);
1495a9f2f41SSatish Balay         ilink = ilink->next;
1501b9fc7fcSBarry Smith       }
151*cf502942SBarry Smith       ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr);
1521b9fc7fcSBarry Smith     }
1531b9fc7fcSBarry Smith   }
1541b9fc7fcSBarry Smith 
15597bbdb24SBarry Smith   if (!jac->pmat) {
156*cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
157*cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
158*cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
159*cf502942SBarry Smith     }
16097bbdb24SBarry Smith   } else {
161*cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
162*cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
163*cf502942SBarry Smith     }
16497bbdb24SBarry Smith   }
16597bbdb24SBarry Smith 
16697bbdb24SBarry Smith   /* set up the individual PCs */
16797bbdb24SBarry Smith   i    = 0;
1685a9f2f41SSatish Balay   ilink = jac->head;
1695a9f2f41SSatish Balay   while (ilink) {
1705a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
1715a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
1725a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
17397bbdb24SBarry Smith     i++;
1745a9f2f41SSatish Balay     ilink = ilink->next;
1750971522cSBarry Smith   }
17697bbdb24SBarry Smith 
17797bbdb24SBarry Smith   /* create work vectors for each split */
1781b9fc7fcSBarry Smith   if (!jac->x) {
17979416396SBarry Smith     Vec xtmp;
18097bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
1815a9f2f41SSatish Balay     ilink = jac->head;
18297bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
18397bbdb24SBarry Smith       Mat A;
1845a9f2f41SSatish Balay       ierr      = KSPGetOperators(ilink->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
1855a9f2f41SSatish Balay       ierr      = MatGetVecs(A,&ilink->x,&ilink->y);CHKERRQ(ierr);
1865a9f2f41SSatish Balay       jac->x[i] = ilink->x;
1875a9f2f41SSatish Balay       jac->y[i] = ilink->y;
1885a9f2f41SSatish Balay       ilink      = ilink->next;
18997bbdb24SBarry Smith     }
19079416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
1911b9fc7fcSBarry Smith 
1925a9f2f41SSatish Balay     ilink = jac->head;
1931b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
1941b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1955a9f2f41SSatish Balay       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
1965a9f2f41SSatish Balay       ilink = ilink->next;
1971b9fc7fcSBarry Smith     }
1981b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
1991b9fc7fcSBarry Smith   }
2000971522cSBarry Smith   PetscFunctionReturn(0);
2010971522cSBarry Smith }
2020971522cSBarry Smith 
2035a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
2045a9f2f41SSatish Balay     (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2055a9f2f41SSatish Balay      VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2065a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
2075a9f2f41SSatish Balay      VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
2085a9f2f41SSatish Balay      VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
20979416396SBarry Smith 
2100971522cSBarry Smith #undef __FUNCT__
2110971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2120971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2130971522cSBarry Smith {
2140971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2150971522cSBarry Smith   PetscErrorCode    ierr;
2165a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
217e25d487eSBarry Smith   PetscInt          bs;
2180971522cSBarry Smith 
2190971522cSBarry Smith   PetscFunctionBegin;
220e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
221e25d487eSBarry Smith   if (bs != jac->bs) {
222e25d487eSBarry Smith     SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector blocksize %D does not match Fieldsplit blocksize %D",bs,jac->bs);
223e25d487eSBarry Smith   }
22479416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2251b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2260971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2275a9f2f41SSatish Balay       while (ilink) {
2285a9f2f41SSatish Balay 	ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2295a9f2f41SSatish Balay 	ilink = ilink->next;
2300971522cSBarry Smith       }
2310971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
2321b9fc7fcSBarry Smith     } else {
2331b9fc7fcSBarry Smith       PetscInt    i = 0;
2341b9fc7fcSBarry Smith 
235efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
2365a9f2f41SSatish Balay       while (ilink) {
2375a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2385a9f2f41SSatish Balay 	ilink = ilink->next;
2391b9fc7fcSBarry Smith 	i++;
2401b9fc7fcSBarry Smith       }
2411b9fc7fcSBarry Smith     }
24279416396SBarry Smith   } else {
24379416396SBarry Smith     if (!jac->w1) {
24479416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
24579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
24679416396SBarry Smith     }
247efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
2485a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2495a9f2f41SSatish Balay     while (ilink->next) {
2505a9f2f41SSatish Balay       ilink = ilink->next;
25179416396SBarry Smith       ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
252efb30889SBarry Smith       ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
2535a9f2f41SSatish Balay       ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
25479416396SBarry Smith     }
25579416396SBarry Smith   }
2560971522cSBarry Smith   PetscFunctionReturn(0);
2570971522cSBarry Smith }
2580971522cSBarry Smith 
2590971522cSBarry Smith #undef __FUNCT__
2600971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
2610971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
2620971522cSBarry Smith {
2630971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2640971522cSBarry Smith   PetscErrorCode    ierr;
2655a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
2660971522cSBarry Smith 
2670971522cSBarry Smith   PetscFunctionBegin;
2685a9f2f41SSatish Balay   while (ilink) {
2695a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
2705a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
2715a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
2725a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
2735a9f2f41SSatish Balay     next = ilink->next;
2745a9f2f41SSatish Balay     ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr);
2755a9f2f41SSatish Balay     ilink = next;
2760971522cSBarry Smith   }
27705b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
27897bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
27997bbdb24SBarry Smith   if (jac->is) {
28097bbdb24SBarry Smith     PetscInt i;
28197bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
28297bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
28397bbdb24SBarry Smith   }
284*cf502942SBarry Smith   if (jac->cis) {
285*cf502942SBarry Smith     PetscInt i;
286*cf502942SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);}
287*cf502942SBarry Smith     ierr = PetscFree(jac->cis);CHKERRQ(ierr);
288*cf502942SBarry Smith   }
28979416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
29079416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
291*cf502942SBarry Smith   ierr = PetscFree(jac->csize);CHKERRQ(ierr);
2920971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
2930971522cSBarry Smith   PetscFunctionReturn(0);
2940971522cSBarry Smith }
2950971522cSBarry Smith 
2960971522cSBarry Smith #undef __FUNCT__
2970971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
2980971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
2990971522cSBarry Smith {
3001b9fc7fcSBarry Smith   PetscErrorCode ierr;
3019dcbbd2bSBarry Smith   PetscInt       i = 0,nfields,fields[12];
3021b9fc7fcSBarry Smith   PetscTruth     flg;
3031b9fc7fcSBarry Smith   char           optionname[128];
3049dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
3051b9fc7fcSBarry Smith 
3060971522cSBarry Smith   PetscFunctionBegin;
3071b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
3089dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
3091b9fc7fcSBarry Smith   while (PETSC_TRUE) {
31013f74950SBarry Smith     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
3111b9fc7fcSBarry Smith     nfields = 12;
3121b9fc7fcSBarry Smith     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
3131b9fc7fcSBarry Smith     if (!flg) break;
3141b9fc7fcSBarry Smith     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
3151b9fc7fcSBarry Smith     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
3161b9fc7fcSBarry Smith   }
3171b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3180971522cSBarry Smith   PetscFunctionReturn(0);
3190971522cSBarry Smith }
3200971522cSBarry Smith 
3210971522cSBarry Smith /*------------------------------------------------------------------------------------*/
3220971522cSBarry Smith 
3230971522cSBarry Smith EXTERN_C_BEGIN
3240971522cSBarry Smith #undef __FUNCT__
3250971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
326dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
3270971522cSBarry Smith {
32897bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3290971522cSBarry Smith   PetscErrorCode    ierr;
3305a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
33169a612a9SBarry Smith   char              prefix[128];
3320971522cSBarry Smith 
3330971522cSBarry Smith   PetscFunctionBegin;
3340971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
3355a9f2f41SSatish Balay   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr);
3365a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
3375a9f2f41SSatish Balay   ilink->nfields = n;
3385a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
3395a9f2f41SSatish Balay   ierr          = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr);
3405a9f2f41SSatish Balay   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
34169a612a9SBarry Smith 
34269a612a9SBarry Smith   if (pc->prefix) {
34313f74950SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
34469a612a9SBarry Smith   } else {
34513f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
34669a612a9SBarry Smith   }
3475a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
3480971522cSBarry Smith 
3490971522cSBarry Smith   if (!next) {
3505a9f2f41SSatish Balay     jac->head = ilink;
3510971522cSBarry Smith   } else {
3520971522cSBarry Smith     while (next->next) {
3530971522cSBarry Smith       next = next->next;
3540971522cSBarry Smith     }
3555a9f2f41SSatish Balay     next->next = ilink;
3560971522cSBarry Smith   }
3570971522cSBarry Smith   jac->nsplits++;
3580971522cSBarry Smith   PetscFunctionReturn(0);
3590971522cSBarry Smith }
3600971522cSBarry Smith EXTERN_C_END
3610971522cSBarry Smith 
3620971522cSBarry Smith 
3630971522cSBarry Smith EXTERN_C_BEGIN
3640971522cSBarry Smith #undef __FUNCT__
36569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
366dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
3670971522cSBarry Smith {
3680971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3690971522cSBarry Smith   PetscErrorCode    ierr;
3700971522cSBarry Smith   PetscInt          cnt = 0;
3715a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
3720971522cSBarry Smith 
3730971522cSBarry Smith   PetscFunctionBegin;
37469a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
3755a9f2f41SSatish Balay   while (ilink) {
3765a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
3775a9f2f41SSatish Balay     ilink = ilink->next;
3780971522cSBarry Smith   }
3790971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
3800971522cSBarry Smith   *n = jac->nsplits;
3810971522cSBarry Smith   PetscFunctionReturn(0);
3820971522cSBarry Smith }
3830971522cSBarry Smith EXTERN_C_END
3840971522cSBarry Smith 
3850971522cSBarry Smith #undef __FUNCT__
3860971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
3870971522cSBarry Smith /*@
3880971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
3890971522cSBarry Smith 
3900971522cSBarry Smith     Collective on PC
3910971522cSBarry Smith 
3920971522cSBarry Smith     Input Parameters:
3930971522cSBarry Smith +   pc  - the preconditioner context
3940971522cSBarry Smith .   n - the number of fields in this split
3950971522cSBarry Smith .   fields - the fields in this split
3960971522cSBarry Smith 
3970971522cSBarry Smith     Level: intermediate
3980971522cSBarry Smith 
39969a612a9SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
4000971522cSBarry Smith 
4010971522cSBarry Smith @*/
402dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
4030971522cSBarry Smith {
4040971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
4050971522cSBarry Smith 
4060971522cSBarry Smith   PetscFunctionBegin;
4070971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4080971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
4090971522cSBarry Smith   if (f) {
4100971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
4110971522cSBarry Smith   }
4120971522cSBarry Smith   PetscFunctionReturn(0);
4130971522cSBarry Smith }
4140971522cSBarry Smith 
4150971522cSBarry Smith #undef __FUNCT__
41669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
4170971522cSBarry Smith /*@C
41869a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
4190971522cSBarry Smith 
42069a612a9SBarry Smith    Collective on KSP
4210971522cSBarry Smith 
4220971522cSBarry Smith    Input Parameter:
4230971522cSBarry Smith .  pc - the preconditioner context
4240971522cSBarry Smith 
4250971522cSBarry Smith    Output Parameters:
4260971522cSBarry Smith +  n - the number of split
42769a612a9SBarry Smith -  pc - the array of KSP contexts
4280971522cSBarry Smith 
4290971522cSBarry Smith    Note:
43069a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
4310971522cSBarry Smith 
43269a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
4330971522cSBarry Smith 
4340971522cSBarry Smith    Level: advanced
4350971522cSBarry Smith 
4360971522cSBarry Smith .seealso: PCFIELDSPLIT
4370971522cSBarry Smith @*/
438dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
4390971522cSBarry Smith {
44069a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
4410971522cSBarry Smith 
4420971522cSBarry Smith   PetscFunctionBegin;
4430971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4440971522cSBarry Smith   PetscValidIntPointer(n,2);
44569a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
4460971522cSBarry Smith   if (f) {
44769a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
4480971522cSBarry Smith   } else {
44969a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
4500971522cSBarry Smith   }
4510971522cSBarry Smith   PetscFunctionReturn(0);
4520971522cSBarry Smith }
4530971522cSBarry Smith 
45479416396SBarry Smith EXTERN_C_BEGIN
45579416396SBarry Smith #undef __FUNCT__
45679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
457dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
45879416396SBarry Smith {
45979416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
46079416396SBarry Smith 
46179416396SBarry Smith   PetscFunctionBegin;
46279416396SBarry Smith   jac->type = type;
46379416396SBarry Smith   PetscFunctionReturn(0);
46479416396SBarry Smith }
46579416396SBarry Smith EXTERN_C_END
46679416396SBarry Smith 
46779416396SBarry Smith #undef __FUNCT__
46879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
46979416396SBarry Smith /*@C
47079416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
47179416396SBarry Smith 
47279416396SBarry Smith    Collective on PC
47379416396SBarry Smith 
47479416396SBarry Smith    Input Parameter:
47579416396SBarry Smith .  pc - the preconditioner context
47679416396SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
47779416396SBarry Smith 
47879416396SBarry Smith    Options Database Key:
47979416396SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
48079416396SBarry Smith 
48179416396SBarry Smith    Level: Developer
48279416396SBarry Smith 
48379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
48479416396SBarry Smith 
48579416396SBarry Smith .seealso: PCCompositeSetType()
48679416396SBarry Smith 
48779416396SBarry Smith @*/
488dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
48979416396SBarry Smith {
49079416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
49179416396SBarry Smith 
49279416396SBarry Smith   PetscFunctionBegin;
49379416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
49479416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
49579416396SBarry Smith   if (f) {
49679416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
49779416396SBarry Smith   }
49879416396SBarry Smith   PetscFunctionReturn(0);
49979416396SBarry Smith }
50079416396SBarry Smith 
5010971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
5020971522cSBarry Smith /*MC
503a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
5040971522cSBarry Smith                   fields or groups of fields
5050971522cSBarry Smith 
5060971522cSBarry Smith 
5070971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
508d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
5090971522cSBarry Smith 
510a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
51169a612a9SBarry Smith          and set the options directly on the resulting KSP object
5120971522cSBarry Smith 
5130971522cSBarry Smith    Level: intermediate
5140971522cSBarry Smith 
51579416396SBarry Smith    Options Database Keys:
5162e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
5172e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
5182e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
5192e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
52079416396SBarry Smith 
5210971522cSBarry Smith    Concepts: physics based preconditioners
5220971522cSBarry Smith 
5230971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
524c8d321e0SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
5250971522cSBarry Smith M*/
5260971522cSBarry Smith 
5270971522cSBarry Smith EXTERN_C_BEGIN
5280971522cSBarry Smith #undef __FUNCT__
5290971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
530dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
5310971522cSBarry Smith {
5320971522cSBarry Smith   PetscErrorCode ierr;
5330971522cSBarry Smith   PC_FieldSplit  *jac;
5340971522cSBarry Smith 
5350971522cSBarry Smith   PetscFunctionBegin;
5360971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
53752e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr);
5380971522cSBarry Smith   jac->bs      = -1;
5390971522cSBarry Smith   jac->nsplits = 0;
54079416396SBarry Smith   jac->type    = PC_COMPOSITE_ADDITIVE;
5410971522cSBarry Smith   pc->data     = (void*)jac;
5420971522cSBarry Smith 
5430971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
5440971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
5450971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
5460971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
5470971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
5480971522cSBarry Smith   pc->ops->applyrichardson   = 0;
5490971522cSBarry Smith 
55069a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
55169a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
5520971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
5530971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
55479416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
55579416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
5560971522cSBarry Smith   PetscFunctionReturn(0);
5570971522cSBarry Smith }
5580971522cSBarry Smith EXTERN_C_END
5590971522cSBarry Smith 
5600971522cSBarry Smith 
561