xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 421e10b8f79bbed49dcc4e803c884835d979c6ea)
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;
1551f519a2SBarry Smith   PC_FieldSplitLink next,previous;
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) {
4251f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);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 
12351f519a2SBarry Smith     bs     = jac->bs;
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         jac->csize[i] = ccsize/nsplit;
13497bbdb24SBarry Smith       } else {
1355a9f2f41SSatish Balay         PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1361b9fc7fcSBarry Smith         PetscTruth sorted;
1375a9f2f41SSatish Balay         ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1381b9fc7fcSBarry Smith         for (j=0; j<nslots; j++) {
1391b9fc7fcSBarry Smith           for (k=0; k<nfields; k++) {
1401b9fc7fcSBarry Smith             ii[nfields*j + k] = rstart + bs*j + fields[k];
14197bbdb24SBarry Smith           }
14297bbdb24SBarry Smith         }
1431b9fc7fcSBarry Smith 	ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr);
144cf502942SBarry Smith         jac->csize[i] = (ccsize/bs)*ilink->nfields;
1451b9fc7fcSBarry Smith         ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr);
1461b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1471b9fc7fcSBarry Smith         ierr = PetscFree(ii);CHKERRQ(ierr);
1485a9f2f41SSatish Balay         ilink = ilink->next;
1491b9fc7fcSBarry Smith       }
150cf502942SBarry Smith       ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr);
1511b9fc7fcSBarry Smith     }
1521b9fc7fcSBarry Smith   }
1531b9fc7fcSBarry Smith 
15497bbdb24SBarry Smith   if (!jac->pmat) {
155cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
156cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
157cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
158cf502942SBarry Smith     }
15997bbdb24SBarry Smith   } else {
160cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
161cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
162cf502942SBarry Smith     }
16397bbdb24SBarry Smith   }
16497bbdb24SBarry Smith 
16597bbdb24SBarry Smith   /* set up the individual PCs */
16697bbdb24SBarry Smith   i    = 0;
1675a9f2f41SSatish Balay   ilink = jac->head;
1685a9f2f41SSatish Balay   while (ilink) {
1695a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
1705a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
1715a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
17297bbdb24SBarry Smith     i++;
1735a9f2f41SSatish Balay     ilink = ilink->next;
1740971522cSBarry Smith   }
17597bbdb24SBarry Smith 
17697bbdb24SBarry Smith   /* create work vectors for each split */
1771b9fc7fcSBarry Smith   if (!jac->x) {
17879416396SBarry Smith     Vec xtmp;
17997bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
1805a9f2f41SSatish Balay     ilink = jac->head;
18197bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
182906ed7ccSBarry Smith       Vec *vl,*vr;
183906ed7ccSBarry Smith 
184906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
185906ed7ccSBarry Smith       ilink->x  = *vr;
186906ed7ccSBarry Smith       ilink->y  = *vl;
187906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
188906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
1895a9f2f41SSatish Balay       jac->x[i] = ilink->x;
1905a9f2f41SSatish Balay       jac->y[i] = ilink->y;
1915a9f2f41SSatish Balay       ilink     = ilink->next;
19297bbdb24SBarry Smith     }
19379416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
1941b9fc7fcSBarry Smith 
1955a9f2f41SSatish Balay     ilink = jac->head;
1961b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
1971b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1985a9f2f41SSatish Balay       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
1995a9f2f41SSatish Balay       ilink = ilink->next;
2001b9fc7fcSBarry Smith     }
2011b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2021b9fc7fcSBarry Smith   }
2030971522cSBarry Smith   PetscFunctionReturn(0);
2040971522cSBarry Smith }
2050971522cSBarry Smith 
2065a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
2075a9f2f41SSatish Balay     (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2085a9f2f41SSatish Balay      VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2095a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
2105a9f2f41SSatish Balay      VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
2115a9f2f41SSatish Balay      VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
21279416396SBarry Smith 
2130971522cSBarry Smith #undef __FUNCT__
2140971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2150971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2160971522cSBarry Smith {
2170971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2180971522cSBarry Smith   PetscErrorCode    ierr;
2195a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
220e25d487eSBarry Smith   PetscInt          bs;
2210971522cSBarry Smith 
2220971522cSBarry Smith   PetscFunctionBegin;
22351f519a2SBarry Smith   CHKMEMQ;
224e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
22551f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
22651f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
22751f519a2SBarry 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 {
237efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
2385a9f2f41SSatish Balay       while (ilink) {
2395a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2405a9f2f41SSatish Balay 	ilink = ilink->next;
2411b9fc7fcSBarry Smith       }
2421b9fc7fcSBarry Smith     }
24379416396SBarry Smith   } else {
24479416396SBarry Smith     if (!jac->w1) {
24579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
24679416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
24779416396SBarry Smith     }
248efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
2495a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2505a9f2f41SSatish Balay     while (ilink->next) {
2515a9f2f41SSatish Balay       ilink = ilink->next;
25279416396SBarry Smith       ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
253efb30889SBarry Smith       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
2545a9f2f41SSatish Balay       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
25579416396SBarry Smith     }
25651f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
25751f519a2SBarry Smith       while (ilink->previous) {
25851f519a2SBarry Smith         ilink = ilink->previous;
25951f519a2SBarry Smith         ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
26051f519a2SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
26151f519a2SBarry Smith         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
26279416396SBarry Smith       }
26351f519a2SBarry Smith     }
26451f519a2SBarry Smith   }
26551f519a2SBarry Smith   CHKMEMQ;
2660971522cSBarry Smith   PetscFunctionReturn(0);
2670971522cSBarry Smith }
2680971522cSBarry Smith 
269*421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
270*421e10b8SBarry Smith     (VecScatterBegin(xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
271*421e10b8SBarry Smith      VecScatterEnd(xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
272*421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
273*421e10b8SBarry Smith      VecScatterBegin(ilink->x,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
274*421e10b8SBarry Smith      VecScatterEnd(ilink->x,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
275*421e10b8SBarry Smith 
276*421e10b8SBarry Smith #undef __FUNCT__
277*421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
278*421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
279*421e10b8SBarry Smith {
280*421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
281*421e10b8SBarry Smith   PetscErrorCode    ierr;
282*421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
283*421e10b8SBarry Smith   PetscInt          bs;
284*421e10b8SBarry Smith 
285*421e10b8SBarry Smith   PetscFunctionBegin;
286*421e10b8SBarry Smith   CHKMEMQ;
287*421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
288*421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
289*421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
290*421e10b8SBarry Smith 
291*421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
292*421e10b8SBarry Smith     if (jac->defaultsplit) {
293*421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
294*421e10b8SBarry Smith       while (ilink) {
295*421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
296*421e10b8SBarry Smith 	ilink = ilink->next;
297*421e10b8SBarry Smith       }
298*421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
299*421e10b8SBarry Smith     } else {
300*421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
301*421e10b8SBarry Smith       while (ilink) {
302*421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
303*421e10b8SBarry Smith 	ilink = ilink->next;
304*421e10b8SBarry Smith       }
305*421e10b8SBarry Smith     }
306*421e10b8SBarry Smith   } else {
307*421e10b8SBarry Smith     if (!jac->w1) {
308*421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
309*421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
310*421e10b8SBarry Smith     }
311*421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
312*421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
313*421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
314*421e10b8SBarry Smith       while (ilink->next) {
315*421e10b8SBarry Smith         ilink = ilink->next;
316*421e10b8SBarry Smith         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
317*421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
318*421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
319*421e10b8SBarry Smith       }
320*421e10b8SBarry Smith       while (ilink->previous) {
321*421e10b8SBarry Smith         ilink = ilink->previous;
322*421e10b8SBarry Smith         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
323*421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
324*421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
325*421e10b8SBarry Smith       }
326*421e10b8SBarry Smith     } else {
327*421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
328*421e10b8SBarry Smith 	ilink = ilink->next;
329*421e10b8SBarry Smith       }
330*421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
331*421e10b8SBarry Smith       while (ilink->previous) {
332*421e10b8SBarry Smith 	ilink = ilink->previous;
333*421e10b8SBarry Smith 	ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
334*421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
335*421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
336*421e10b8SBarry Smith       }
337*421e10b8SBarry Smith     }
338*421e10b8SBarry Smith   }
339*421e10b8SBarry Smith   CHKMEMQ;
340*421e10b8SBarry Smith   PetscFunctionReturn(0);
341*421e10b8SBarry Smith }
342*421e10b8SBarry Smith 
3430971522cSBarry Smith #undef __FUNCT__
3440971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
3450971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
3460971522cSBarry Smith {
3470971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3480971522cSBarry Smith   PetscErrorCode    ierr;
3495a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
3500971522cSBarry Smith 
3510971522cSBarry Smith   PetscFunctionBegin;
3525a9f2f41SSatish Balay   while (ilink) {
3535a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
3545a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
3555a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
3565a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
3575a9f2f41SSatish Balay     next = ilink->next;
3585a9f2f41SSatish Balay     ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr);
3595a9f2f41SSatish Balay     ilink = next;
3600971522cSBarry Smith   }
36105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
36297bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
36397bbdb24SBarry Smith   if (jac->is) {
36497bbdb24SBarry Smith     PetscInt i;
36597bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
36697bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
36797bbdb24SBarry Smith   }
368cf502942SBarry Smith   if (jac->cis) {
369cf502942SBarry Smith     PetscInt i;
370cf502942SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);}
371cf502942SBarry Smith     ierr = PetscFree(jac->cis);CHKERRQ(ierr);
372cf502942SBarry Smith   }
37379416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
37479416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
375cf502942SBarry Smith   ierr = PetscFree(jac->csize);CHKERRQ(ierr);
3760971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
3770971522cSBarry Smith   PetscFunctionReturn(0);
3780971522cSBarry Smith }
3790971522cSBarry Smith 
3800971522cSBarry Smith #undef __FUNCT__
3810971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
3820971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
3830971522cSBarry Smith {
3841b9fc7fcSBarry Smith   PetscErrorCode ierr;
38551f519a2SBarry Smith   PetscInt       i = 0,nfields,*fields,bs;
3861b9fc7fcSBarry Smith   PetscTruth     flg;
3871b9fc7fcSBarry Smith   char           optionname[128];
3889dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
3891b9fc7fcSBarry Smith 
3900971522cSBarry Smith   PetscFunctionBegin;
3911b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
39251f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
39351f519a2SBarry Smith   if (flg) {
39451f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
39551f519a2SBarry Smith   }
39651f519a2SBarry Smith   if (jac->bs <= 0) {
39751f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr);
39851f519a2SBarry Smith   }
3999dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
40051f519a2SBarry Smith   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
4011b9fc7fcSBarry Smith   while (PETSC_TRUE) {
40213f74950SBarry Smith     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
40351f519a2SBarry Smith     nfields = jac->bs;
4041b9fc7fcSBarry Smith     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
4051b9fc7fcSBarry Smith     if (!flg) break;
4061b9fc7fcSBarry Smith     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
4071b9fc7fcSBarry Smith     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
4081b9fc7fcSBarry Smith   }
40951f519a2SBarry Smith   ierr = PetscFree(fields);CHKERRQ(ierr);
4101b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4110971522cSBarry Smith   PetscFunctionReturn(0);
4120971522cSBarry Smith }
4130971522cSBarry Smith 
4140971522cSBarry Smith /*------------------------------------------------------------------------------------*/
4150971522cSBarry Smith 
4160971522cSBarry Smith EXTERN_C_BEGIN
4170971522cSBarry Smith #undef __FUNCT__
4180971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
419dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
4200971522cSBarry Smith {
42197bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4220971522cSBarry Smith   PetscErrorCode    ierr;
4235a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
42469a612a9SBarry Smith   char              prefix[128];
42551f519a2SBarry Smith   PetscInt          i;
4260971522cSBarry Smith 
4270971522cSBarry Smith   PetscFunctionBegin;
4280971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
42951f519a2SBarry Smith   for (i=0; i<n; i++) {
43051f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
43151f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
43251f519a2SBarry Smith   }
4335a9f2f41SSatish Balay   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr);
4345a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
4355a9f2f41SSatish Balay   ilink->nfields = n;
4365a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
4375a9f2f41SSatish Balay   ierr           = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr);
4385a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
43969a612a9SBarry Smith 
44069a612a9SBarry Smith   if (pc->prefix) {
44113f74950SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
44269a612a9SBarry Smith   } else {
44313f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
44469a612a9SBarry Smith   }
4455a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
4460971522cSBarry Smith 
4470971522cSBarry Smith   if (!next) {
4485a9f2f41SSatish Balay     jac->head       = ilink;
44951f519a2SBarry Smith     ilink->previous = PETSC_NULL;
4500971522cSBarry Smith   } else {
4510971522cSBarry Smith     while (next->next) {
4520971522cSBarry Smith       next = next->next;
4530971522cSBarry Smith     }
4545a9f2f41SSatish Balay     next->next      = ilink;
45551f519a2SBarry Smith     ilink->previous = next;
4560971522cSBarry Smith   }
4570971522cSBarry Smith   jac->nsplits++;
4580971522cSBarry Smith   PetscFunctionReturn(0);
4590971522cSBarry Smith }
4600971522cSBarry Smith EXTERN_C_END
4610971522cSBarry Smith 
4620971522cSBarry Smith 
4630971522cSBarry Smith EXTERN_C_BEGIN
4640971522cSBarry Smith #undef __FUNCT__
46569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
466dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
4670971522cSBarry Smith {
4680971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4690971522cSBarry Smith   PetscErrorCode    ierr;
4700971522cSBarry Smith   PetscInt          cnt = 0;
4715a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
4720971522cSBarry Smith 
4730971522cSBarry Smith   PetscFunctionBegin;
47469a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
4755a9f2f41SSatish Balay   while (ilink) {
4765a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
4775a9f2f41SSatish Balay     ilink = ilink->next;
4780971522cSBarry Smith   }
4790971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
4800971522cSBarry Smith   *n = jac->nsplits;
4810971522cSBarry Smith   PetscFunctionReturn(0);
4820971522cSBarry Smith }
4830971522cSBarry Smith EXTERN_C_END
4840971522cSBarry Smith 
4850971522cSBarry Smith #undef __FUNCT__
4860971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
4870971522cSBarry Smith /*@
4880971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
4890971522cSBarry Smith 
4900971522cSBarry Smith     Collective on PC
4910971522cSBarry Smith 
4920971522cSBarry Smith     Input Parameters:
4930971522cSBarry Smith +   pc  - the preconditioner context
4940971522cSBarry Smith .   n - the number of fields in this split
4950971522cSBarry Smith .   fields - the fields in this split
4960971522cSBarry Smith 
4970971522cSBarry Smith     Level: intermediate
4980971522cSBarry Smith 
49951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
5000971522cSBarry Smith 
5010971522cSBarry Smith @*/
502dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
5030971522cSBarry Smith {
5040971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
5050971522cSBarry Smith 
5060971522cSBarry Smith   PetscFunctionBegin;
5070971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5080971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
5090971522cSBarry Smith   if (f) {
5100971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
5110971522cSBarry Smith   }
5120971522cSBarry Smith   PetscFunctionReturn(0);
5130971522cSBarry Smith }
5140971522cSBarry Smith 
5150971522cSBarry Smith #undef __FUNCT__
51651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
51751f519a2SBarry Smith /*@
51851f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
51951f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
52051f519a2SBarry Smith 
52151f519a2SBarry Smith     Collective on PC
52251f519a2SBarry Smith 
52351f519a2SBarry Smith     Input Parameters:
52451f519a2SBarry Smith +   pc  - the preconditioner context
52551f519a2SBarry Smith -   bs - the block size
52651f519a2SBarry Smith 
52751f519a2SBarry Smith     Level: intermediate
52851f519a2SBarry Smith 
52951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
53051f519a2SBarry Smith 
53151f519a2SBarry Smith @*/
53251f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
53351f519a2SBarry Smith {
53451f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
53551f519a2SBarry Smith 
53651f519a2SBarry Smith   PetscFunctionBegin;
53751f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
53851f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
53951f519a2SBarry Smith   if (f) {
54051f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
54151f519a2SBarry Smith   }
54251f519a2SBarry Smith   PetscFunctionReturn(0);
54351f519a2SBarry Smith }
54451f519a2SBarry Smith 
54551f519a2SBarry Smith #undef __FUNCT__
54669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
5470971522cSBarry Smith /*@C
54869a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
5490971522cSBarry Smith 
55069a612a9SBarry Smith    Collective on KSP
5510971522cSBarry Smith 
5520971522cSBarry Smith    Input Parameter:
5530971522cSBarry Smith .  pc - the preconditioner context
5540971522cSBarry Smith 
5550971522cSBarry Smith    Output Parameters:
5560971522cSBarry Smith +  n - the number of split
55769a612a9SBarry Smith -  pc - the array of KSP contexts
5580971522cSBarry Smith 
5590971522cSBarry Smith    Note:
56069a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
5610971522cSBarry Smith 
56269a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
5630971522cSBarry Smith 
5640971522cSBarry Smith    Level: advanced
5650971522cSBarry Smith 
5660971522cSBarry Smith .seealso: PCFIELDSPLIT
5670971522cSBarry Smith @*/
568dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
5690971522cSBarry Smith {
57069a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
5710971522cSBarry Smith 
5720971522cSBarry Smith   PetscFunctionBegin;
5730971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5740971522cSBarry Smith   PetscValidIntPointer(n,2);
57569a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
5760971522cSBarry Smith   if (f) {
57769a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
5780971522cSBarry Smith   } else {
57969a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
5800971522cSBarry Smith   }
5810971522cSBarry Smith   PetscFunctionReturn(0);
5820971522cSBarry Smith }
5830971522cSBarry Smith 
58479416396SBarry Smith EXTERN_C_BEGIN
58579416396SBarry Smith #undef __FUNCT__
58679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
587dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
58879416396SBarry Smith {
58979416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
59079416396SBarry Smith 
59179416396SBarry Smith   PetscFunctionBegin;
59279416396SBarry Smith   jac->type = type;
59379416396SBarry Smith   PetscFunctionReturn(0);
59479416396SBarry Smith }
59579416396SBarry Smith EXTERN_C_END
59679416396SBarry Smith 
59751f519a2SBarry Smith EXTERN_C_BEGIN
59851f519a2SBarry Smith #undef __FUNCT__
59951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
60051f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
60151f519a2SBarry Smith {
60251f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
60351f519a2SBarry Smith 
60451f519a2SBarry Smith   PetscFunctionBegin;
60551f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
60651f519a2SBarry Smith   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
60751f519a2SBarry Smith   jac->bs = bs;
60851f519a2SBarry Smith   PetscFunctionReturn(0);
60951f519a2SBarry Smith }
61051f519a2SBarry Smith EXTERN_C_END
61151f519a2SBarry Smith 
61279416396SBarry Smith #undef __FUNCT__
61379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
61479416396SBarry Smith /*@C
61579416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
61679416396SBarry Smith 
61779416396SBarry Smith    Collective on PC
61879416396SBarry Smith 
61979416396SBarry Smith    Input Parameter:
62079416396SBarry Smith .  pc - the preconditioner context
62179416396SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
62279416396SBarry Smith 
62379416396SBarry Smith    Options Database Key:
62479416396SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
62579416396SBarry Smith 
62679416396SBarry Smith    Level: Developer
62779416396SBarry Smith 
62879416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
62979416396SBarry Smith 
63079416396SBarry Smith .seealso: PCCompositeSetType()
63179416396SBarry Smith 
63279416396SBarry Smith @*/
633dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
63479416396SBarry Smith {
63579416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
63679416396SBarry Smith 
63779416396SBarry Smith   PetscFunctionBegin;
63879416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
63979416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
64079416396SBarry Smith   if (f) {
64179416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
64279416396SBarry Smith   }
64379416396SBarry Smith   PetscFunctionReturn(0);
64479416396SBarry Smith }
64579416396SBarry Smith 
6460971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
6470971522cSBarry Smith /*MC
648a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
6490971522cSBarry Smith                   fields or groups of fields
6500971522cSBarry Smith 
6510971522cSBarry Smith 
6520971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
653d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
6540971522cSBarry Smith 
655a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
65669a612a9SBarry Smith          and set the options directly on the resulting KSP object
6570971522cSBarry Smith 
6580971522cSBarry Smith    Level: intermediate
6590971522cSBarry Smith 
66079416396SBarry Smith    Options Database Keys:
6612e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
6622e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
6632e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
66451f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
6652e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
66679416396SBarry Smith 
6670971522cSBarry Smith    Concepts: physics based preconditioners
6680971522cSBarry Smith 
6690971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
670c8d321e0SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
6710971522cSBarry Smith M*/
6720971522cSBarry Smith 
6730971522cSBarry Smith EXTERN_C_BEGIN
6740971522cSBarry Smith #undef __FUNCT__
6750971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
676dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
6770971522cSBarry Smith {
6780971522cSBarry Smith   PetscErrorCode ierr;
6790971522cSBarry Smith   PC_FieldSplit  *jac;
6800971522cSBarry Smith 
6810971522cSBarry Smith   PetscFunctionBegin;
6820971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
68352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr);
6840971522cSBarry Smith   jac->bs        = -1;
6850971522cSBarry Smith   jac->nsplits   = 0;
68679416396SBarry Smith   jac->type      = PC_COMPOSITE_ADDITIVE;
68751f519a2SBarry Smith 
6880971522cSBarry Smith   pc->data     = (void*)jac;
6890971522cSBarry Smith 
6900971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
691*421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
6920971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
6930971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
6940971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
6950971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
6960971522cSBarry Smith   pc->ops->applyrichardson   = 0;
6970971522cSBarry Smith 
69869a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
69969a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
7000971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
7010971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
70279416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
70379416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
70451f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
70551f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
7060971522cSBarry Smith   PetscFunctionReturn(0);
7070971522cSBarry Smith }
7080971522cSBarry Smith EXTERN_C_END
7090971522cSBarry Smith 
7100971522cSBarry Smith 
711