xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 906ed7cc33fecbafab01746eec64dcdcc8a4842f)
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;
22cf502942SBarry Smith   PetscInt          nsplits,*csize;
2379416396SBarry Smith   Vec               *x,*y,w1,w2;
2497bbdb24SBarry Smith   Mat               *pmat;
25cf502942SBarry 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) {
42cf502942SBarry 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;
111cf502942SBarry 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);
125cf502942SBarry 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);
128cf502942SBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->cis);CHKERRQ(ierr);
129cf502942SBarry 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);
133cf502942SBarry 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);
134cf502942SBarry 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);
145cf502942SBarry 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       }
151cf502942SBarry Smith       ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr);
1521b9fc7fcSBarry Smith     }
1531b9fc7fcSBarry Smith   }
1541b9fc7fcSBarry Smith 
15597bbdb24SBarry Smith   if (!jac->pmat) {
156cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
157cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
158cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
159cf502942SBarry Smith     }
16097bbdb24SBarry Smith   } else {
161cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
162cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
163cf502942SBarry 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++) {
183*906ed7ccSBarry Smith       Vec *vl,*vr;
184*906ed7ccSBarry Smith 
185*906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
186*906ed7ccSBarry Smith       ilink->x  = *vr;
187*906ed7ccSBarry Smith       ilink->y  = *vl;
188*906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
189*906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
1905a9f2f41SSatish Balay       jac->x[i] = ilink->x;
1915a9f2f41SSatish Balay       jac->y[i] = ilink->y;
1925a9f2f41SSatish Balay       ilink     = ilink->next;
19397bbdb24SBarry Smith     }
19479416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
1951b9fc7fcSBarry Smith 
1965a9f2f41SSatish Balay     ilink = jac->head;
1971b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
1981b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1995a9f2f41SSatish Balay       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
2005a9f2f41SSatish Balay       ilink = ilink->next;
2011b9fc7fcSBarry Smith     }
2021b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2031b9fc7fcSBarry Smith   }
2040971522cSBarry Smith   PetscFunctionReturn(0);
2050971522cSBarry Smith }
2060971522cSBarry Smith 
2075a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
2085a9f2f41SSatish Balay     (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2095a9f2f41SSatish Balay      VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2105a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
2115a9f2f41SSatish Balay      VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
2125a9f2f41SSatish Balay      VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
21379416396SBarry Smith 
2140971522cSBarry Smith #undef __FUNCT__
2150971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2160971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2170971522cSBarry Smith {
2180971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2190971522cSBarry Smith   PetscErrorCode    ierr;
2205a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
221e25d487eSBarry Smith   PetscInt          bs;
2220971522cSBarry Smith 
2230971522cSBarry Smith   PetscFunctionBegin;
224e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
225e25d487eSBarry Smith   if (bs != jac->bs) {
226e25d487eSBarry Smith     SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector blocksize %D does not match Fieldsplit blocksize %D",bs,jac->bs);
227e25d487eSBarry Smith   }
22879416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2291b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2300971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2315a9f2f41SSatish Balay       while (ilink) {
2325a9f2f41SSatish Balay 	ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2335a9f2f41SSatish Balay 	ilink = ilink->next;
2340971522cSBarry Smith       }
2350971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
2361b9fc7fcSBarry Smith     } else {
2371b9fc7fcSBarry Smith       PetscInt    i = 0;
2381b9fc7fcSBarry Smith 
239efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
2405a9f2f41SSatish Balay       while (ilink) {
2415a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2425a9f2f41SSatish Balay 	ilink = ilink->next;
2431b9fc7fcSBarry Smith 	i++;
2441b9fc7fcSBarry Smith       }
2451b9fc7fcSBarry Smith     }
24679416396SBarry Smith   } else {
24779416396SBarry Smith     if (!jac->w1) {
24879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
24979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
25079416396SBarry Smith     }
251efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
2525a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2535a9f2f41SSatish Balay     while (ilink->next) {
2545a9f2f41SSatish Balay       ilink = ilink->next;
25579416396SBarry Smith       ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
256efb30889SBarry Smith       ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
2575a9f2f41SSatish Balay       ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
25879416396SBarry Smith     }
25979416396SBarry Smith   }
2600971522cSBarry Smith   PetscFunctionReturn(0);
2610971522cSBarry Smith }
2620971522cSBarry Smith 
2630971522cSBarry Smith #undef __FUNCT__
2640971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
2650971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
2660971522cSBarry Smith {
2670971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2680971522cSBarry Smith   PetscErrorCode    ierr;
2695a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
2700971522cSBarry Smith 
2710971522cSBarry Smith   PetscFunctionBegin;
2725a9f2f41SSatish Balay   while (ilink) {
2735a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
2745a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
2755a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
2765a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
2775a9f2f41SSatish Balay     next = ilink->next;
2785a9f2f41SSatish Balay     ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr);
2795a9f2f41SSatish Balay     ilink = next;
2800971522cSBarry Smith   }
28105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
28297bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
28397bbdb24SBarry Smith   if (jac->is) {
28497bbdb24SBarry Smith     PetscInt i;
28597bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
28697bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
28797bbdb24SBarry Smith   }
288cf502942SBarry Smith   if (jac->cis) {
289cf502942SBarry Smith     PetscInt i;
290cf502942SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);}
291cf502942SBarry Smith     ierr = PetscFree(jac->cis);CHKERRQ(ierr);
292cf502942SBarry Smith   }
29379416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
29479416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
295cf502942SBarry Smith   ierr = PetscFree(jac->csize);CHKERRQ(ierr);
2960971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
2970971522cSBarry Smith   PetscFunctionReturn(0);
2980971522cSBarry Smith }
2990971522cSBarry Smith 
3000971522cSBarry Smith #undef __FUNCT__
3010971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
3020971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
3030971522cSBarry Smith {
3041b9fc7fcSBarry Smith   PetscErrorCode ierr;
3059dcbbd2bSBarry Smith   PetscInt       i = 0,nfields,fields[12];
3061b9fc7fcSBarry Smith   PetscTruth     flg;
3071b9fc7fcSBarry Smith   char           optionname[128];
3089dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
3091b9fc7fcSBarry Smith 
3100971522cSBarry Smith   PetscFunctionBegin;
3111b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
3129dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
3131b9fc7fcSBarry Smith   while (PETSC_TRUE) {
31413f74950SBarry Smith     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
3151b9fc7fcSBarry Smith     nfields = 12;
3161b9fc7fcSBarry Smith     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
3171b9fc7fcSBarry Smith     if (!flg) break;
3181b9fc7fcSBarry Smith     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
3191b9fc7fcSBarry Smith     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
3201b9fc7fcSBarry Smith   }
3211b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3220971522cSBarry Smith   PetscFunctionReturn(0);
3230971522cSBarry Smith }
3240971522cSBarry Smith 
3250971522cSBarry Smith /*------------------------------------------------------------------------------------*/
3260971522cSBarry Smith 
3270971522cSBarry Smith EXTERN_C_BEGIN
3280971522cSBarry Smith #undef __FUNCT__
3290971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
330dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
3310971522cSBarry Smith {
33297bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3330971522cSBarry Smith   PetscErrorCode    ierr;
3345a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
33569a612a9SBarry Smith   char              prefix[128];
3360971522cSBarry Smith 
3370971522cSBarry Smith   PetscFunctionBegin;
3380971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
3395a9f2f41SSatish Balay   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr);
3405a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
3415a9f2f41SSatish Balay   ilink->nfields = n;
3425a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
3435a9f2f41SSatish Balay   ierr          = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr);
3445a9f2f41SSatish Balay   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
34569a612a9SBarry Smith 
34669a612a9SBarry Smith   if (pc->prefix) {
34713f74950SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
34869a612a9SBarry Smith   } else {
34913f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
35069a612a9SBarry Smith   }
3515a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
3520971522cSBarry Smith 
3530971522cSBarry Smith   if (!next) {
3545a9f2f41SSatish Balay     jac->head = ilink;
3550971522cSBarry Smith   } else {
3560971522cSBarry Smith     while (next->next) {
3570971522cSBarry Smith       next = next->next;
3580971522cSBarry Smith     }
3595a9f2f41SSatish Balay     next->next = ilink;
3600971522cSBarry Smith   }
3610971522cSBarry Smith   jac->nsplits++;
3620971522cSBarry Smith   PetscFunctionReturn(0);
3630971522cSBarry Smith }
3640971522cSBarry Smith EXTERN_C_END
3650971522cSBarry Smith 
3660971522cSBarry Smith 
3670971522cSBarry Smith EXTERN_C_BEGIN
3680971522cSBarry Smith #undef __FUNCT__
36969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
370dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
3710971522cSBarry Smith {
3720971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3730971522cSBarry Smith   PetscErrorCode    ierr;
3740971522cSBarry Smith   PetscInt          cnt = 0;
3755a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
3760971522cSBarry Smith 
3770971522cSBarry Smith   PetscFunctionBegin;
37869a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
3795a9f2f41SSatish Balay   while (ilink) {
3805a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
3815a9f2f41SSatish Balay     ilink = ilink->next;
3820971522cSBarry Smith   }
3830971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
3840971522cSBarry Smith   *n = jac->nsplits;
3850971522cSBarry Smith   PetscFunctionReturn(0);
3860971522cSBarry Smith }
3870971522cSBarry Smith EXTERN_C_END
3880971522cSBarry Smith 
3890971522cSBarry Smith #undef __FUNCT__
3900971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
3910971522cSBarry Smith /*@
3920971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
3930971522cSBarry Smith 
3940971522cSBarry Smith     Collective on PC
3950971522cSBarry Smith 
3960971522cSBarry Smith     Input Parameters:
3970971522cSBarry Smith +   pc  - the preconditioner context
3980971522cSBarry Smith .   n - the number of fields in this split
3990971522cSBarry Smith .   fields - the fields in this split
4000971522cSBarry Smith 
4010971522cSBarry Smith     Level: intermediate
4020971522cSBarry Smith 
40369a612a9SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
4040971522cSBarry Smith 
4050971522cSBarry Smith @*/
406dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
4070971522cSBarry Smith {
4080971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
4090971522cSBarry Smith 
4100971522cSBarry Smith   PetscFunctionBegin;
4110971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4120971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
4130971522cSBarry Smith   if (f) {
4140971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
4150971522cSBarry Smith   }
4160971522cSBarry Smith   PetscFunctionReturn(0);
4170971522cSBarry Smith }
4180971522cSBarry Smith 
4190971522cSBarry Smith #undef __FUNCT__
42069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
4210971522cSBarry Smith /*@C
42269a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
4230971522cSBarry Smith 
42469a612a9SBarry Smith    Collective on KSP
4250971522cSBarry Smith 
4260971522cSBarry Smith    Input Parameter:
4270971522cSBarry Smith .  pc - the preconditioner context
4280971522cSBarry Smith 
4290971522cSBarry Smith    Output Parameters:
4300971522cSBarry Smith +  n - the number of split
43169a612a9SBarry Smith -  pc - the array of KSP contexts
4320971522cSBarry Smith 
4330971522cSBarry Smith    Note:
43469a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
4350971522cSBarry Smith 
43669a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
4370971522cSBarry Smith 
4380971522cSBarry Smith    Level: advanced
4390971522cSBarry Smith 
4400971522cSBarry Smith .seealso: PCFIELDSPLIT
4410971522cSBarry Smith @*/
442dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
4430971522cSBarry Smith {
44469a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
4450971522cSBarry Smith 
4460971522cSBarry Smith   PetscFunctionBegin;
4470971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4480971522cSBarry Smith   PetscValidIntPointer(n,2);
44969a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
4500971522cSBarry Smith   if (f) {
45169a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
4520971522cSBarry Smith   } else {
45369a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
4540971522cSBarry Smith   }
4550971522cSBarry Smith   PetscFunctionReturn(0);
4560971522cSBarry Smith }
4570971522cSBarry Smith 
45879416396SBarry Smith EXTERN_C_BEGIN
45979416396SBarry Smith #undef __FUNCT__
46079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
461dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
46279416396SBarry Smith {
46379416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
46479416396SBarry Smith 
46579416396SBarry Smith   PetscFunctionBegin;
46679416396SBarry Smith   jac->type = type;
46779416396SBarry Smith   PetscFunctionReturn(0);
46879416396SBarry Smith }
46979416396SBarry Smith EXTERN_C_END
47079416396SBarry Smith 
47179416396SBarry Smith #undef __FUNCT__
47279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
47379416396SBarry Smith /*@C
47479416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
47579416396SBarry Smith 
47679416396SBarry Smith    Collective on PC
47779416396SBarry Smith 
47879416396SBarry Smith    Input Parameter:
47979416396SBarry Smith .  pc - the preconditioner context
48079416396SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
48179416396SBarry Smith 
48279416396SBarry Smith    Options Database Key:
48379416396SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
48479416396SBarry Smith 
48579416396SBarry Smith    Level: Developer
48679416396SBarry Smith 
48779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
48879416396SBarry Smith 
48979416396SBarry Smith .seealso: PCCompositeSetType()
49079416396SBarry Smith 
49179416396SBarry Smith @*/
492dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
49379416396SBarry Smith {
49479416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
49579416396SBarry Smith 
49679416396SBarry Smith   PetscFunctionBegin;
49779416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
49879416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
49979416396SBarry Smith   if (f) {
50079416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
50179416396SBarry Smith   }
50279416396SBarry Smith   PetscFunctionReturn(0);
50379416396SBarry Smith }
50479416396SBarry Smith 
5050971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
5060971522cSBarry Smith /*MC
507a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
5080971522cSBarry Smith                   fields or groups of fields
5090971522cSBarry Smith 
5100971522cSBarry Smith 
5110971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
512d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
5130971522cSBarry Smith 
514a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
51569a612a9SBarry Smith          and set the options directly on the resulting KSP object
5160971522cSBarry Smith 
5170971522cSBarry Smith    Level: intermediate
5180971522cSBarry Smith 
51979416396SBarry Smith    Options Database Keys:
5202e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
5212e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
5222e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
5232e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
52479416396SBarry Smith 
5250971522cSBarry Smith    Concepts: physics based preconditioners
5260971522cSBarry Smith 
5270971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
528c8d321e0SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
5290971522cSBarry Smith M*/
5300971522cSBarry Smith 
5310971522cSBarry Smith EXTERN_C_BEGIN
5320971522cSBarry Smith #undef __FUNCT__
5330971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
534dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
5350971522cSBarry Smith {
5360971522cSBarry Smith   PetscErrorCode ierr;
5370971522cSBarry Smith   PC_FieldSplit  *jac;
5380971522cSBarry Smith 
5390971522cSBarry Smith   PetscFunctionBegin;
5400971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
54152e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr);
5420971522cSBarry Smith   jac->bs      = -1;
5430971522cSBarry Smith   jac->nsplits = 0;
54479416396SBarry Smith   jac->type    = PC_COMPOSITE_ADDITIVE;
5450971522cSBarry Smith   pc->data     = (void*)jac;
5460971522cSBarry Smith 
5470971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
5480971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
5490971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
5500971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
5510971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
5520971522cSBarry Smith   pc->ops->applyrichardson   = 0;
5530971522cSBarry Smith 
55469a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
55569a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
5560971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
5570971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
55879416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
55979416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
5600971522cSBarry Smith   PetscFunctionReturn(0);
5610971522cSBarry Smith }
5620971522cSBarry Smith EXTERN_C_END
5630971522cSBarry Smith 
5640971522cSBarry Smith 
565