xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 704ba839cd484dc2900d48fe4907f851bfd76db1)
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;
15*704ba839SBarry Smith   IS                is,cis;
16*704ba839SBarry Smith   PetscInt          csize;
1751f519a2SBarry Smith   PC_FieldSplitLink next,previous;
180971522cSBarry Smith };
190971522cSBarry Smith 
200971522cSBarry Smith typedef struct {
2179416396SBarry Smith   PCCompositeType   type;              /* additive or multiplicative */
2297bbdb24SBarry Smith   PetscTruth        defaultsplit;
230971522cSBarry Smith   PetscInt          bs;
24*704ba839SBarry Smith   PetscInt          nsplits;
2579416396SBarry Smith   Vec               *x,*y,w1,w2;
2697bbdb24SBarry Smith   Mat               *pmat;
27*704ba839SBarry Smith   PetscTruth        issetup;
2897bbdb24SBarry Smith   PC_FieldSplitLink head;
290971522cSBarry Smith } PC_FieldSplit;
300971522cSBarry Smith 
310971522cSBarry Smith #undef __FUNCT__
320971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
330971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
340971522cSBarry Smith {
350971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
360971522cSBarry Smith   PetscErrorCode    ierr;
370971522cSBarry Smith   PetscTruth        iascii;
380971522cSBarry Smith   PetscInt          i,j;
395a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
400971522cSBarry Smith 
410971522cSBarry Smith   PetscFunctionBegin;
420971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
430971522cSBarry Smith   if (iascii) {
4451f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
4569a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
460971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
470971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
480971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
4979416396SBarry Smith       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
505a9f2f41SSatish Balay       for (j=0; j<ilink->nfields; j++) {
5179416396SBarry Smith         if (j > 0) {
5279416396SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
5379416396SBarry Smith         }
545a9f2f41SSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
550971522cSBarry Smith       }
560971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
5779416396SBarry Smith       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
585a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
595a9f2f41SSatish Balay       ilink = ilink->next;
600971522cSBarry Smith     }
610971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
620971522cSBarry Smith   } else {
630971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
640971522cSBarry Smith   }
650971522cSBarry Smith   PetscFunctionReturn(0);
660971522cSBarry Smith }
670971522cSBarry Smith 
680971522cSBarry Smith #undef __FUNCT__
6969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
7069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
710971522cSBarry Smith {
720971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
730971522cSBarry Smith   PetscErrorCode    ierr;
745a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
7569a612a9SBarry Smith   PetscInt          i;
7679416396SBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields;
770971522cSBarry Smith 
780971522cSBarry Smith   PetscFunctionBegin;
797adad957SLisandro Dalcin   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
805a9f2f41SSatish Balay   if (!ilink || flg) {
81ae15b995SBarry Smith     ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
82521d7252SBarry Smith     if (jac->bs <= 0) {
83*704ba839SBarry Smith       if (pc->pmat) {
84521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
85*704ba839SBarry Smith       } else {
86*704ba839SBarry Smith         jac->bs = 1;
87*704ba839SBarry Smith       }
88521d7252SBarry Smith     }
8979416396SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
9079416396SBarry Smith     ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
915a9f2f41SSatish Balay     while (ilink) {
925a9f2f41SSatish Balay       for (i=0; i<ilink->nfields; i++) {
935a9f2f41SSatish Balay         fields[ilink->fields[i]] = PETSC_TRUE;
94521d7252SBarry Smith       }
955a9f2f41SSatish Balay       ilink = ilink->next;
9679416396SBarry Smith     }
9797bbdb24SBarry Smith     jac->defaultsplit = PETSC_TRUE;
9879416396SBarry Smith     for (i=0; i<jac->bs; i++) {
9979416396SBarry Smith       if (!fields[i]) {
10079416396SBarry Smith 	ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
10179416396SBarry Smith       } else {
10279416396SBarry Smith         jac->defaultsplit = PETSC_FALSE;
10379416396SBarry Smith       }
10479416396SBarry Smith     }
105521d7252SBarry Smith   }
10669a612a9SBarry Smith   PetscFunctionReturn(0);
10769a612a9SBarry Smith }
10869a612a9SBarry Smith 
10969a612a9SBarry Smith 
11069a612a9SBarry Smith #undef __FUNCT__
11169a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
11269a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
11369a612a9SBarry Smith {
11469a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11569a612a9SBarry Smith   PetscErrorCode    ierr;
1165a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
117cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
11869a612a9SBarry Smith   MatStructure      flag = pc->flag;
119ccb205f8SBarry Smith   PetscTruth        sorted;
12069a612a9SBarry Smith 
12169a612a9SBarry Smith   PetscFunctionBegin;
12269a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
12397bbdb24SBarry Smith   nsplit = jac->nsplits;
1245a9f2f41SSatish Balay   ilink  = jac->head;
12597bbdb24SBarry Smith 
12697bbdb24SBarry Smith   /* get the matrices for each split */
127*704ba839SBarry Smith   if (!jac->issetup) {
1281b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
12997bbdb24SBarry Smith 
130*704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
131*704ba839SBarry Smith 
132*704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
13351f519a2SBarry Smith     bs     = jac->bs;
13497bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
135cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
1361b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1371b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1381b9fc7fcSBarry Smith       if (jac->defaultsplit) {
139*704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
140*704ba839SBarry Smith         ilink->csize = ccsize/nsplit;
141*704ba839SBarry Smith       } else if (!ilink->is) {
142ccb205f8SBarry Smith         if (ilink->nfields > 1) {
1435a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1445a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1451b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
1461b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
1471b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
14897bbdb24SBarry Smith             }
14997bbdb24SBarry Smith           }
150*704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
151ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
152ccb205f8SBarry Smith         } else {
153*704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
154ccb205f8SBarry Smith         }
155*704ba839SBarry Smith         ilink->csize = (ccsize/bs)*ilink->nfields;
156*704ba839SBarry Smith         ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
1571b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1581b9fc7fcSBarry Smith       }
159*704ba839SBarry Smith       ierr  = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
160*704ba839SBarry Smith       ilink = ilink->next;
1611b9fc7fcSBarry Smith     }
1621b9fc7fcSBarry Smith   }
1631b9fc7fcSBarry Smith 
164*704ba839SBarry Smith   ilink  = jac->head;
16597bbdb24SBarry Smith   if (!jac->pmat) {
166cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
167cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
168*704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
169*704ba839SBarry Smith       ilink = ilink->next;
170cf502942SBarry Smith     }
17197bbdb24SBarry Smith   } else {
172cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
173*704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
174*704ba839SBarry Smith       ilink = ilink->next;
175cf502942SBarry Smith     }
17697bbdb24SBarry Smith   }
17797bbdb24SBarry Smith 
17897bbdb24SBarry Smith   /* set up the individual PCs */
17997bbdb24SBarry Smith   i    = 0;
1805a9f2f41SSatish Balay   ilink = jac->head;
1815a9f2f41SSatish Balay   while (ilink) {
1825a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
1835a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
1845a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
18597bbdb24SBarry Smith     i++;
1865a9f2f41SSatish Balay     ilink = ilink->next;
1870971522cSBarry Smith   }
18897bbdb24SBarry Smith 
18997bbdb24SBarry Smith   /* create work vectors for each split */
1901b9fc7fcSBarry Smith   if (!jac->x) {
19179416396SBarry Smith     Vec xtmp;
19297bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
1935a9f2f41SSatish Balay     ilink = jac->head;
19497bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
195906ed7ccSBarry Smith       Vec *vl,*vr;
196906ed7ccSBarry Smith 
197906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
198906ed7ccSBarry Smith       ilink->x  = *vr;
199906ed7ccSBarry Smith       ilink->y  = *vl;
200906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
201906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
2025a9f2f41SSatish Balay       jac->x[i] = ilink->x;
2035a9f2f41SSatish Balay       jac->y[i] = ilink->y;
2045a9f2f41SSatish Balay       ilink     = ilink->next;
20597bbdb24SBarry Smith     }
20679416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
2071b9fc7fcSBarry Smith 
2085a9f2f41SSatish Balay     ilink = jac->head;
2091b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
2101b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
211*704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
2125a9f2f41SSatish Balay       ilink = ilink->next;
2131b9fc7fcSBarry Smith     }
2141b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2151b9fc7fcSBarry Smith   }
2160971522cSBarry Smith   PetscFunctionReturn(0);
2170971522cSBarry Smith }
2180971522cSBarry Smith 
2195a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
220ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
221ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
2225a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
223ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
224ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
22579416396SBarry Smith 
2260971522cSBarry Smith #undef __FUNCT__
2270971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2280971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2290971522cSBarry Smith {
2300971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2310971522cSBarry Smith   PetscErrorCode    ierr;
2325a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
233e25d487eSBarry Smith   PetscInt          bs;
2340971522cSBarry Smith 
2350971522cSBarry Smith   PetscFunctionBegin;
23651f519a2SBarry Smith   CHKMEMQ;
237e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
23851f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
23951f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
24051f519a2SBarry Smith 
24179416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2421b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2430971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2445a9f2f41SSatish Balay       while (ilink) {
2455a9f2f41SSatish Balay 	ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2465a9f2f41SSatish Balay 	ilink = ilink->next;
2470971522cSBarry Smith       }
2480971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
2491b9fc7fcSBarry Smith     } else {
250efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
2515a9f2f41SSatish Balay       while (ilink) {
2525a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2535a9f2f41SSatish Balay 	ilink = ilink->next;
2541b9fc7fcSBarry Smith       }
2551b9fc7fcSBarry Smith     }
25679416396SBarry Smith   } else {
25779416396SBarry Smith     if (!jac->w1) {
25879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
25979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
26079416396SBarry Smith     }
261efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
2625a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2635a9f2f41SSatish Balay     while (ilink->next) {
2645a9f2f41SSatish Balay       ilink = ilink->next;
26579416396SBarry Smith       ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
266efb30889SBarry Smith       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
2675a9f2f41SSatish Balay       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
26879416396SBarry Smith     }
26951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
27051f519a2SBarry Smith       while (ilink->previous) {
27151f519a2SBarry Smith         ilink = ilink->previous;
27251f519a2SBarry Smith         ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
27351f519a2SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
27451f519a2SBarry Smith         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
27579416396SBarry Smith       }
27651f519a2SBarry Smith     }
27751f519a2SBarry Smith   }
27851f519a2SBarry Smith   CHKMEMQ;
2790971522cSBarry Smith   PetscFunctionReturn(0);
2800971522cSBarry Smith }
2810971522cSBarry Smith 
282421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
283ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
284ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
285421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
286ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
287ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
288421e10b8SBarry Smith 
289421e10b8SBarry Smith #undef __FUNCT__
290421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
291421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
292421e10b8SBarry Smith {
293421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
294421e10b8SBarry Smith   PetscErrorCode    ierr;
295421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
296421e10b8SBarry Smith   PetscInt          bs;
297421e10b8SBarry Smith 
298421e10b8SBarry Smith   PetscFunctionBegin;
299421e10b8SBarry Smith   CHKMEMQ;
300421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
301421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
302421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
303421e10b8SBarry Smith 
304421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
305421e10b8SBarry Smith     if (jac->defaultsplit) {
306421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
307421e10b8SBarry Smith       while (ilink) {
308421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
309421e10b8SBarry Smith 	ilink = ilink->next;
310421e10b8SBarry Smith       }
311421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
312421e10b8SBarry Smith     } else {
313421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
314421e10b8SBarry Smith       while (ilink) {
315421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
316421e10b8SBarry Smith 	ilink = ilink->next;
317421e10b8SBarry Smith       }
318421e10b8SBarry Smith     }
319421e10b8SBarry Smith   } else {
320421e10b8SBarry Smith     if (!jac->w1) {
321421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
322421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
323421e10b8SBarry Smith     }
324421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
325421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
326421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
327421e10b8SBarry Smith       while (ilink->next) {
328421e10b8SBarry Smith         ilink = ilink->next;
329421e10b8SBarry Smith         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
330421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
331421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
332421e10b8SBarry Smith       }
333421e10b8SBarry Smith       while (ilink->previous) {
334421e10b8SBarry Smith         ilink = ilink->previous;
335421e10b8SBarry Smith         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
336421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
337421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
338421e10b8SBarry Smith       }
339421e10b8SBarry Smith     } else {
340421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
341421e10b8SBarry Smith 	ilink = ilink->next;
342421e10b8SBarry Smith       }
343421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
344421e10b8SBarry Smith       while (ilink->previous) {
345421e10b8SBarry Smith 	ilink = ilink->previous;
346421e10b8SBarry Smith 	ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
347421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
348421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
349421e10b8SBarry Smith       }
350421e10b8SBarry Smith     }
351421e10b8SBarry Smith   }
352421e10b8SBarry Smith   CHKMEMQ;
353421e10b8SBarry Smith   PetscFunctionReturn(0);
354421e10b8SBarry Smith }
355421e10b8SBarry Smith 
3560971522cSBarry Smith #undef __FUNCT__
3570971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
3580971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
3590971522cSBarry Smith {
3600971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3610971522cSBarry Smith   PetscErrorCode    ierr;
3625a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
3630971522cSBarry Smith 
3640971522cSBarry Smith   PetscFunctionBegin;
3655a9f2f41SSatish Balay   while (ilink) {
3665a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
3675a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
3685a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
3695a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
370*704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
371*704ba839SBarry Smith     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
3725a9f2f41SSatish Balay     next = ilink->next;
373*704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
374*704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
3755a9f2f41SSatish Balay     ilink = next;
3760971522cSBarry Smith   }
37705b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
37897bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
37979416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
38079416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
3810971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
3820971522cSBarry Smith   PetscFunctionReturn(0);
3830971522cSBarry Smith }
3840971522cSBarry Smith 
3850971522cSBarry Smith #undef __FUNCT__
3860971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
3870971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
3880971522cSBarry Smith {
3891b9fc7fcSBarry Smith   PetscErrorCode ierr;
39051f519a2SBarry Smith   PetscInt       i = 0,nfields,*fields,bs;
3911b9fc7fcSBarry Smith   PetscTruth     flg;
3921b9fc7fcSBarry Smith   char           optionname[128];
3939dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
3941b9fc7fcSBarry Smith 
3950971522cSBarry Smith   PetscFunctionBegin;
3961b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
39751f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
39851f519a2SBarry Smith   if (flg) {
39951f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
40051f519a2SBarry Smith   }
40151f519a2SBarry Smith   if (jac->bs <= 0) {
402*704ba839SBarry Smith 
40351f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr);
40451f519a2SBarry Smith   }
4059dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
40651f519a2SBarry Smith   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
4071b9fc7fcSBarry Smith   while (PETSC_TRUE) {
40813f74950SBarry Smith     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
40951f519a2SBarry Smith     nfields = jac->bs;
4101b9fc7fcSBarry Smith     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
4111b9fc7fcSBarry Smith     if (!flg) break;
4121b9fc7fcSBarry Smith     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
4131b9fc7fcSBarry Smith     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
4141b9fc7fcSBarry Smith   }
41551f519a2SBarry Smith   ierr = PetscFree(fields);CHKERRQ(ierr);
4161b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4170971522cSBarry Smith   PetscFunctionReturn(0);
4180971522cSBarry Smith }
4190971522cSBarry Smith 
4200971522cSBarry Smith /*------------------------------------------------------------------------------------*/
4210971522cSBarry Smith 
4220971522cSBarry Smith EXTERN_C_BEGIN
4230971522cSBarry Smith #undef __FUNCT__
4240971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
425dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
4260971522cSBarry Smith {
42797bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4280971522cSBarry Smith   PetscErrorCode    ierr;
4295a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
43069a612a9SBarry Smith   char              prefix[128];
43151f519a2SBarry Smith   PetscInt          i;
4320971522cSBarry Smith 
4330971522cSBarry Smith   PetscFunctionBegin;
4340971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
43551f519a2SBarry Smith   for (i=0; i<n; i++) {
43651f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
43751f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
43851f519a2SBarry Smith   }
439*704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
440*704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
4415a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
4425a9f2f41SSatish Balay   ilink->nfields = n;
4435a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
4447adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
4455a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
44669a612a9SBarry Smith 
4477adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
4487adad957SLisandro Dalcin     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
44969a612a9SBarry Smith   } else {
45013f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
45169a612a9SBarry Smith   }
4525a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
4530971522cSBarry Smith 
4540971522cSBarry Smith   if (!next) {
4555a9f2f41SSatish Balay     jac->head       = ilink;
45651f519a2SBarry Smith     ilink->previous = PETSC_NULL;
4570971522cSBarry Smith   } else {
4580971522cSBarry Smith     while (next->next) {
4590971522cSBarry Smith       next = next->next;
4600971522cSBarry Smith     }
4615a9f2f41SSatish Balay     next->next      = ilink;
46251f519a2SBarry Smith     ilink->previous = next;
4630971522cSBarry Smith   }
4640971522cSBarry Smith   jac->nsplits++;
4650971522cSBarry Smith   PetscFunctionReturn(0);
4660971522cSBarry Smith }
4670971522cSBarry Smith EXTERN_C_END
4680971522cSBarry Smith 
4690971522cSBarry Smith 
4700971522cSBarry Smith EXTERN_C_BEGIN
4710971522cSBarry Smith #undef __FUNCT__
47269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
473dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
4740971522cSBarry Smith {
4750971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4760971522cSBarry Smith   PetscErrorCode    ierr;
4770971522cSBarry Smith   PetscInt          cnt = 0;
4785a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
4790971522cSBarry Smith 
4800971522cSBarry Smith   PetscFunctionBegin;
48169a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
4825a9f2f41SSatish Balay   while (ilink) {
4835a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
4845a9f2f41SSatish Balay     ilink = ilink->next;
4850971522cSBarry Smith   }
4860971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
4870971522cSBarry Smith   *n = jac->nsplits;
4880971522cSBarry Smith   PetscFunctionReturn(0);
4890971522cSBarry Smith }
4900971522cSBarry Smith EXTERN_C_END
4910971522cSBarry Smith 
492*704ba839SBarry Smith EXTERN_C_BEGIN
493*704ba839SBarry Smith #undef __FUNCT__
494*704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
495*704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
496*704ba839SBarry Smith {
497*704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
498*704ba839SBarry Smith   PetscErrorCode    ierr;
499*704ba839SBarry Smith     PC_FieldSplitLink ilink, next = jac->head;
500*704ba839SBarry Smith   char              prefix[128];
501*704ba839SBarry Smith 
502*704ba839SBarry Smith   PetscFunctionBegin;
503*704ba839SBarry Smith   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,0,PetscInt,&ilink->fields);CHKERRQ(ierr);
504*704ba839SBarry Smith   ilink->nfields = -1;
505*704ba839SBarry Smith   ilink->next    = PETSC_NULL;
506*704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
507*704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
508*704ba839SBarry Smith 
509*704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
510*704ba839SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
511*704ba839SBarry Smith   } else {
512*704ba839SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
513*704ba839SBarry Smith   }
514*704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
515*704ba839SBarry Smith 
516*704ba839SBarry Smith   if (!next) {
517*704ba839SBarry Smith     jac->head       = ilink;
518*704ba839SBarry Smith     ilink->previous = PETSC_NULL;
519*704ba839SBarry Smith   } else {
520*704ba839SBarry Smith     while (next->next) {
521*704ba839SBarry Smith       next = next->next;
522*704ba839SBarry Smith     }
523*704ba839SBarry Smith     next->next      = ilink;
524*704ba839SBarry Smith     ilink->previous = next;
525*704ba839SBarry Smith   }
526*704ba839SBarry Smith   jac->nsplits++;
527*704ba839SBarry Smith 
528*704ba839SBarry Smith   PetscFunctionReturn(0);
529*704ba839SBarry Smith }
530*704ba839SBarry Smith EXTERN_C_END
531*704ba839SBarry Smith 
5320971522cSBarry Smith #undef __FUNCT__
5330971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
5340971522cSBarry Smith /*@
5350971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
5360971522cSBarry Smith 
5370971522cSBarry Smith     Collective on PC
5380971522cSBarry Smith 
5390971522cSBarry Smith     Input Parameters:
5400971522cSBarry Smith +   pc  - the preconditioner context
5410971522cSBarry Smith .   n - the number of fields in this split
5420971522cSBarry Smith .   fields - the fields in this split
5430971522cSBarry Smith 
5440971522cSBarry Smith     Level: intermediate
5450971522cSBarry Smith 
54651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
5470971522cSBarry Smith 
5480971522cSBarry Smith @*/
549dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
5500971522cSBarry Smith {
5510971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
5520971522cSBarry Smith 
5530971522cSBarry Smith   PetscFunctionBegin;
5540971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5550971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
5560971522cSBarry Smith   if (f) {
5570971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
5580971522cSBarry Smith   }
5590971522cSBarry Smith   PetscFunctionReturn(0);
5600971522cSBarry Smith }
5610971522cSBarry Smith 
5620971522cSBarry Smith #undef __FUNCT__
563*704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
564*704ba839SBarry Smith /*@
565*704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
566*704ba839SBarry Smith 
567*704ba839SBarry Smith     Collective on PC
568*704ba839SBarry Smith 
569*704ba839SBarry Smith     Input Parameters:
570*704ba839SBarry Smith +   pc  - the preconditioner context
571*704ba839SBarry Smith .   is - the index set that defines the vector elements in this field
572*704ba839SBarry Smith 
573*704ba839SBarry Smith     Level: intermediate
574*704ba839SBarry Smith 
575*704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
576*704ba839SBarry Smith 
577*704ba839SBarry Smith @*/
578*704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
579*704ba839SBarry Smith {
580*704ba839SBarry Smith   PetscErrorCode ierr,(*f)(PC,IS);
581*704ba839SBarry Smith 
582*704ba839SBarry Smith   PetscFunctionBegin;
583*704ba839SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
584*704ba839SBarry Smith   PetscValidHeaderSpecific(is,IS_COOKIE,1);
585*704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
586*704ba839SBarry Smith   if (f) {
587*704ba839SBarry Smith     ierr = (*f)(pc,is);CHKERRQ(ierr);
588*704ba839SBarry Smith   }
589*704ba839SBarry Smith   PetscFunctionReturn(0);
590*704ba839SBarry Smith }
591*704ba839SBarry Smith 
592*704ba839SBarry Smith #undef __FUNCT__
59351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
59451f519a2SBarry Smith /*@
59551f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
59651f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
59751f519a2SBarry Smith 
59851f519a2SBarry Smith     Collective on PC
59951f519a2SBarry Smith 
60051f519a2SBarry Smith     Input Parameters:
60151f519a2SBarry Smith +   pc  - the preconditioner context
60251f519a2SBarry Smith -   bs - the block size
60351f519a2SBarry Smith 
60451f519a2SBarry Smith     Level: intermediate
60551f519a2SBarry Smith 
60651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
60751f519a2SBarry Smith 
60851f519a2SBarry Smith @*/
60951f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
61051f519a2SBarry Smith {
61151f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
61251f519a2SBarry Smith 
61351f519a2SBarry Smith   PetscFunctionBegin;
61451f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
61551f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
61651f519a2SBarry Smith   if (f) {
61751f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
61851f519a2SBarry Smith   }
61951f519a2SBarry Smith   PetscFunctionReturn(0);
62051f519a2SBarry Smith }
62151f519a2SBarry Smith 
62251f519a2SBarry Smith #undef __FUNCT__
62369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
6240971522cSBarry Smith /*@C
62569a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
6260971522cSBarry Smith 
62769a612a9SBarry Smith    Collective on KSP
6280971522cSBarry Smith 
6290971522cSBarry Smith    Input Parameter:
6300971522cSBarry Smith .  pc - the preconditioner context
6310971522cSBarry Smith 
6320971522cSBarry Smith    Output Parameters:
6330971522cSBarry Smith +  n - the number of split
63469a612a9SBarry Smith -  pc - the array of KSP contexts
6350971522cSBarry Smith 
6360971522cSBarry Smith    Note:
63769a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
6380971522cSBarry Smith 
63969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
6400971522cSBarry Smith 
6410971522cSBarry Smith    Level: advanced
6420971522cSBarry Smith 
6430971522cSBarry Smith .seealso: PCFIELDSPLIT
6440971522cSBarry Smith @*/
645dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
6460971522cSBarry Smith {
64769a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
6480971522cSBarry Smith 
6490971522cSBarry Smith   PetscFunctionBegin;
6500971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6510971522cSBarry Smith   PetscValidIntPointer(n,2);
65269a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
6530971522cSBarry Smith   if (f) {
65469a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
6550971522cSBarry Smith   } else {
65669a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
6570971522cSBarry Smith   }
6580971522cSBarry Smith   PetscFunctionReturn(0);
6590971522cSBarry Smith }
6600971522cSBarry Smith 
66179416396SBarry Smith EXTERN_C_BEGIN
66279416396SBarry Smith #undef __FUNCT__
66379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
664dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
66579416396SBarry Smith {
66679416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
66779416396SBarry Smith 
66879416396SBarry Smith   PetscFunctionBegin;
66979416396SBarry Smith   jac->type = type;
67079416396SBarry Smith   PetscFunctionReturn(0);
67179416396SBarry Smith }
67279416396SBarry Smith EXTERN_C_END
67379416396SBarry Smith 
67451f519a2SBarry Smith EXTERN_C_BEGIN
67551f519a2SBarry Smith #undef __FUNCT__
67651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
67751f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
67851f519a2SBarry Smith {
67951f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
68051f519a2SBarry Smith 
68151f519a2SBarry Smith   PetscFunctionBegin;
68251f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
68351f519a2SBarry 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);
68451f519a2SBarry Smith   jac->bs = bs;
68551f519a2SBarry Smith   PetscFunctionReturn(0);
68651f519a2SBarry Smith }
68751f519a2SBarry Smith EXTERN_C_END
68851f519a2SBarry Smith 
68979416396SBarry Smith #undef __FUNCT__
69079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
691bc08b0f1SBarry Smith /*@
69279416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
69379416396SBarry Smith 
69479416396SBarry Smith    Collective on PC
69579416396SBarry Smith 
69679416396SBarry Smith    Input Parameter:
69779416396SBarry Smith .  pc - the preconditioner context
698ce9499c7SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
69979416396SBarry Smith 
70079416396SBarry Smith    Options Database Key:
701ce9499c7SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
70279416396SBarry Smith 
70379416396SBarry Smith    Level: Developer
70479416396SBarry Smith 
70579416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
70679416396SBarry Smith 
70779416396SBarry Smith .seealso: PCCompositeSetType()
70879416396SBarry Smith 
70979416396SBarry Smith @*/
710dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
71179416396SBarry Smith {
71279416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
71379416396SBarry Smith 
71479416396SBarry Smith   PetscFunctionBegin;
71579416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
71679416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
71779416396SBarry Smith   if (f) {
71879416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
71979416396SBarry Smith   }
72079416396SBarry Smith   PetscFunctionReturn(0);
72179416396SBarry Smith }
72279416396SBarry Smith 
7230971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
7240971522cSBarry Smith /*MC
725a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
7260971522cSBarry Smith                   fields or groups of fields
7270971522cSBarry Smith 
7280971522cSBarry Smith 
7290971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
730d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
7310971522cSBarry Smith 
732a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
73369a612a9SBarry Smith          and set the options directly on the resulting KSP object
7340971522cSBarry Smith 
7350971522cSBarry Smith    Level: intermediate
7360971522cSBarry Smith 
73779416396SBarry Smith    Options Database Keys:
7382e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
7392e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
7402e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
74151f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
7422e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
74379416396SBarry Smith 
7440971522cSBarry Smith    Concepts: physics based preconditioners
7450971522cSBarry Smith 
7460971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
747c8d321e0SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
7480971522cSBarry Smith M*/
7490971522cSBarry Smith 
7500971522cSBarry Smith EXTERN_C_BEGIN
7510971522cSBarry Smith #undef __FUNCT__
7520971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
753dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
7540971522cSBarry Smith {
7550971522cSBarry Smith   PetscErrorCode ierr;
7560971522cSBarry Smith   PC_FieldSplit  *jac;
7570971522cSBarry Smith 
7580971522cSBarry Smith   PetscFunctionBegin;
75938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
7600971522cSBarry Smith   jac->bs        = -1;
7610971522cSBarry Smith   jac->nsplits   = 0;
76279416396SBarry Smith   jac->type      = PC_COMPOSITE_ADDITIVE;
76351f519a2SBarry Smith 
7640971522cSBarry Smith   pc->data     = (void*)jac;
7650971522cSBarry Smith 
7660971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
767421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
7680971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
7690971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
7700971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
7710971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
7720971522cSBarry Smith   pc->ops->applyrichardson   = 0;
7730971522cSBarry Smith 
77469a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
77569a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
7760971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
7770971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
778*704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
779*704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
78079416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
78179416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
78251f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
78351f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
7840971522cSBarry Smith   PetscFunctionReturn(0);
7850971522cSBarry Smith }
7860971522cSBarry Smith EXTERN_C_END
7870971522cSBarry Smith 
7880971522cSBarry Smith 
789