xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision ccb205f84456404d0833d80f1cbd41e068105bb1)
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;
113*ccb205f8SBarry Smith   PetscTruth        sorted;
11469a612a9SBarry Smith 
11569a612a9SBarry Smith   PetscFunctionBegin;
11669a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
11797bbdb24SBarry Smith   nsplit = jac->nsplits;
1185a9f2f41SSatish Balay   ilink  = jac->head;
11997bbdb24SBarry Smith 
12097bbdb24SBarry Smith   /* get the matrices for each split */
12197bbdb24SBarry Smith   if (!jac->is) {
1221b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
12397bbdb24SBarry Smith 
12451f519a2SBarry Smith     bs     = jac->bs;
12597bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
126cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
1271b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1281b9fc7fcSBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr);
129cf502942SBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->cis);CHKERRQ(ierr);
130cf502942SBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(PetscInt),&jac->csize);CHKERRQ(ierr);
1311b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1321b9fc7fcSBarry Smith       if (jac->defaultsplit) {
1331b9fc7fcSBarry Smith 	ierr     = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr);
134cf502942SBarry Smith         jac->csize[i] = ccsize/nsplit;
13597bbdb24SBarry Smith       } else {
136*ccb205f8SBarry Smith 	printf("ilink nfields %d\n",ilink->nfields);
137*ccb205f8SBarry Smith         if (ilink->nfields > 1) {
1385a9f2f41SSatish Balay 	  PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1395a9f2f41SSatish Balay 	  ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1401b9fc7fcSBarry Smith 	  for (j=0; j<nslots; j++) {
1411b9fc7fcSBarry Smith 	    for (k=0; k<nfields; k++) {
1421b9fc7fcSBarry Smith 	      ii[nfields*j + k] = rstart + bs*j + fields[k];
14397bbdb24SBarry Smith 	    }
14497bbdb24SBarry Smith 	  }
1451b9fc7fcSBarry Smith 	  ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr);
146*ccb205f8SBarry Smith 	  ierr = PetscFree(ii);CHKERRQ(ierr);
147*ccb205f8SBarry Smith         } else {
148*ccb205f8SBarry Smith           ierr = ISCreateStride(pc->comm,nslots,ilink->fields[0],bs,&jac->is[i]);CHKERRQ(ierr);
149*ccb205f8SBarry Smith         }
150cf502942SBarry Smith         jac->csize[i] = (ccsize/bs)*ilink->nfields;
1511b9fc7fcSBarry Smith         ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr);
1521b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1535a9f2f41SSatish Balay         ilink = ilink->next;
1541b9fc7fcSBarry Smith       }
155cf502942SBarry Smith       ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr);
1561b9fc7fcSBarry Smith     }
1571b9fc7fcSBarry Smith   }
1581b9fc7fcSBarry Smith 
15997bbdb24SBarry Smith   if (!jac->pmat) {
160cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
161cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
162cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
163cf502942SBarry Smith     }
16497bbdb24SBarry Smith   } else {
165cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
166cf502942SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
167cf502942SBarry Smith     }
16897bbdb24SBarry Smith   }
16997bbdb24SBarry Smith 
17097bbdb24SBarry Smith   /* set up the individual PCs */
17197bbdb24SBarry Smith   i    = 0;
1725a9f2f41SSatish Balay   ilink = jac->head;
1735a9f2f41SSatish Balay   while (ilink) {
1745a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
1755a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
1765a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
17797bbdb24SBarry Smith     i++;
1785a9f2f41SSatish Balay     ilink = ilink->next;
1790971522cSBarry Smith   }
18097bbdb24SBarry Smith 
18197bbdb24SBarry Smith   /* create work vectors for each split */
1821b9fc7fcSBarry Smith   if (!jac->x) {
18379416396SBarry Smith     Vec xtmp;
18497bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
1855a9f2f41SSatish Balay     ilink = jac->head;
18697bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
187906ed7ccSBarry Smith       Vec *vl,*vr;
188906ed7ccSBarry Smith 
189906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
190906ed7ccSBarry Smith       ilink->x  = *vr;
191906ed7ccSBarry Smith       ilink->y  = *vl;
192906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
193906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
1945a9f2f41SSatish Balay       jac->x[i] = ilink->x;
1955a9f2f41SSatish Balay       jac->y[i] = ilink->y;
1965a9f2f41SSatish Balay       ilink     = ilink->next;
19797bbdb24SBarry Smith     }
19879416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
1991b9fc7fcSBarry Smith 
2005a9f2f41SSatish Balay     ilink = jac->head;
2011b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
2021b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
2035a9f2f41SSatish Balay       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
2045a9f2f41SSatish Balay       ilink = ilink->next;
2051b9fc7fcSBarry Smith     }
2061b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2071b9fc7fcSBarry Smith   }
2080971522cSBarry Smith   PetscFunctionReturn(0);
2090971522cSBarry Smith }
2100971522cSBarry Smith 
2115a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
2125a9f2f41SSatish Balay     (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2135a9f2f41SSatish Balay      VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
2145a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
2155a9f2f41SSatish Balay      VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
2165a9f2f41SSatish Balay      VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
21779416396SBarry Smith 
2180971522cSBarry Smith #undef __FUNCT__
2190971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2200971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2210971522cSBarry Smith {
2220971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2230971522cSBarry Smith   PetscErrorCode    ierr;
2245a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
225e25d487eSBarry Smith   PetscInt          bs;
2260971522cSBarry Smith 
2270971522cSBarry Smith   PetscFunctionBegin;
22851f519a2SBarry Smith   CHKMEMQ;
229e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
23051f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
23151f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
23251f519a2SBarry Smith 
23379416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2341b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2350971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2365a9f2f41SSatish Balay       while (ilink) {
2375a9f2f41SSatish Balay 	ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2385a9f2f41SSatish Balay 	ilink = ilink->next;
2390971522cSBarry Smith       }
2400971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
2411b9fc7fcSBarry Smith     } else {
242efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
2435a9f2f41SSatish Balay       while (ilink) {
2445a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2455a9f2f41SSatish Balay 	ilink = ilink->next;
2461b9fc7fcSBarry Smith       }
2471b9fc7fcSBarry Smith     }
24879416396SBarry Smith   } else {
24979416396SBarry Smith     if (!jac->w1) {
25079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
25179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
25279416396SBarry Smith     }
253efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
2545a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2555a9f2f41SSatish Balay     while (ilink->next) {
2565a9f2f41SSatish Balay       ilink = ilink->next;
25779416396SBarry Smith       ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
258efb30889SBarry Smith       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
2595a9f2f41SSatish Balay       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
26079416396SBarry Smith     }
26151f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
26251f519a2SBarry Smith       while (ilink->previous) {
26351f519a2SBarry Smith         ilink = ilink->previous;
26451f519a2SBarry Smith         ierr  = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
26551f519a2SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
26651f519a2SBarry Smith         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
26779416396SBarry Smith       }
26851f519a2SBarry Smith     }
26951f519a2SBarry Smith   }
27051f519a2SBarry Smith   CHKMEMQ;
2710971522cSBarry Smith   PetscFunctionReturn(0);
2720971522cSBarry Smith }
2730971522cSBarry Smith 
274421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
275421e10b8SBarry Smith     (VecScatterBegin(xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
276421e10b8SBarry Smith      VecScatterEnd(xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
277421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
278421e10b8SBarry Smith      VecScatterBegin(ilink->x,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
279421e10b8SBarry Smith      VecScatterEnd(ilink->x,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
280421e10b8SBarry Smith 
281421e10b8SBarry Smith #undef __FUNCT__
282421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
283421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
284421e10b8SBarry Smith {
285421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
286421e10b8SBarry Smith   PetscErrorCode    ierr;
287421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
288421e10b8SBarry Smith   PetscInt          bs;
289421e10b8SBarry Smith 
290421e10b8SBarry Smith   PetscFunctionBegin;
291421e10b8SBarry Smith   CHKMEMQ;
292421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
293421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
294421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
295421e10b8SBarry Smith 
296421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
297421e10b8SBarry Smith     if (jac->defaultsplit) {
298421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
299421e10b8SBarry Smith       while (ilink) {
300421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
301421e10b8SBarry Smith 	ilink = ilink->next;
302421e10b8SBarry Smith       }
303421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
304421e10b8SBarry Smith     } else {
305421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
306421e10b8SBarry Smith       while (ilink) {
307421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
308421e10b8SBarry Smith 	ilink = ilink->next;
309421e10b8SBarry Smith       }
310421e10b8SBarry Smith     }
311421e10b8SBarry Smith   } else {
312421e10b8SBarry Smith     if (!jac->w1) {
313421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
314421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
315421e10b8SBarry Smith     }
316421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
317421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
318421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
319421e10b8SBarry Smith       while (ilink->next) {
320421e10b8SBarry Smith         ilink = ilink->next;
321421e10b8SBarry Smith         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
322421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
323421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
324421e10b8SBarry Smith       }
325421e10b8SBarry Smith       while (ilink->previous) {
326421e10b8SBarry Smith         ilink = ilink->previous;
327421e10b8SBarry Smith         ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
328421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
329421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
330421e10b8SBarry Smith       }
331421e10b8SBarry Smith     } else {
332421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
333421e10b8SBarry Smith 	ilink = ilink->next;
334421e10b8SBarry Smith       }
335421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
336421e10b8SBarry Smith       while (ilink->previous) {
337421e10b8SBarry Smith 	ilink = ilink->previous;
338421e10b8SBarry Smith 	ierr  = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr);
339421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
340421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
341421e10b8SBarry Smith       }
342421e10b8SBarry Smith     }
343421e10b8SBarry Smith   }
344421e10b8SBarry Smith   CHKMEMQ;
345421e10b8SBarry Smith   PetscFunctionReturn(0);
346421e10b8SBarry Smith }
347421e10b8SBarry Smith 
3480971522cSBarry Smith #undef __FUNCT__
3490971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
3500971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
3510971522cSBarry Smith {
3520971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3530971522cSBarry Smith   PetscErrorCode    ierr;
3545a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
3550971522cSBarry Smith 
3560971522cSBarry Smith   PetscFunctionBegin;
3575a9f2f41SSatish Balay   while (ilink) {
3585a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
3595a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
3605a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
3615a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
3625a9f2f41SSatish Balay     next = ilink->next;
3635a9f2f41SSatish Balay     ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr);
3645a9f2f41SSatish Balay     ilink = next;
3650971522cSBarry Smith   }
36605b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
36797bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
36897bbdb24SBarry Smith   if (jac->is) {
36997bbdb24SBarry Smith     PetscInt i;
37097bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
37197bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
37297bbdb24SBarry Smith   }
373cf502942SBarry Smith   if (jac->cis) {
374cf502942SBarry Smith     PetscInt i;
375cf502942SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);}
376cf502942SBarry Smith     ierr = PetscFree(jac->cis);CHKERRQ(ierr);
377cf502942SBarry Smith   }
37879416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
37979416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
380cf502942SBarry Smith   ierr = PetscFree(jac->csize);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) {
40251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr);
40351f519a2SBarry Smith   }
4049dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
40551f519a2SBarry Smith   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
4061b9fc7fcSBarry Smith   while (PETSC_TRUE) {
40713f74950SBarry Smith     sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
40851f519a2SBarry Smith     nfields = jac->bs;
4091b9fc7fcSBarry Smith     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
4101b9fc7fcSBarry Smith     if (!flg) break;
4111b9fc7fcSBarry Smith     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
4121b9fc7fcSBarry Smith     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
4131b9fc7fcSBarry Smith   }
41451f519a2SBarry Smith   ierr = PetscFree(fields);CHKERRQ(ierr);
4151b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4160971522cSBarry Smith   PetscFunctionReturn(0);
4170971522cSBarry Smith }
4180971522cSBarry Smith 
4190971522cSBarry Smith /*------------------------------------------------------------------------------------*/
4200971522cSBarry Smith 
4210971522cSBarry Smith EXTERN_C_BEGIN
4220971522cSBarry Smith #undef __FUNCT__
4230971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
424dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
4250971522cSBarry Smith {
42697bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4270971522cSBarry Smith   PetscErrorCode    ierr;
4285a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
42969a612a9SBarry Smith   char              prefix[128];
43051f519a2SBarry Smith   PetscInt          i;
4310971522cSBarry Smith 
4320971522cSBarry Smith   PetscFunctionBegin;
4330971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
43451f519a2SBarry Smith   for (i=0; i<n; i++) {
43551f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
43651f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
43751f519a2SBarry Smith   }
4385a9f2f41SSatish Balay   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr);
4395a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
4405a9f2f41SSatish Balay   ilink->nfields = n;
4415a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
4425a9f2f41SSatish Balay   ierr           = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr);
4435a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
44469a612a9SBarry Smith 
44569a612a9SBarry Smith   if (pc->prefix) {
44613f74950SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
44769a612a9SBarry Smith   } else {
44813f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
44969a612a9SBarry Smith   }
4505a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
4510971522cSBarry Smith 
4520971522cSBarry Smith   if (!next) {
4535a9f2f41SSatish Balay     jac->head       = ilink;
45451f519a2SBarry Smith     ilink->previous = PETSC_NULL;
4550971522cSBarry Smith   } else {
4560971522cSBarry Smith     while (next->next) {
4570971522cSBarry Smith       next = next->next;
4580971522cSBarry Smith     }
4595a9f2f41SSatish Balay     next->next      = ilink;
46051f519a2SBarry Smith     ilink->previous = next;
4610971522cSBarry Smith   }
4620971522cSBarry Smith   jac->nsplits++;
4630971522cSBarry Smith   PetscFunctionReturn(0);
4640971522cSBarry Smith }
4650971522cSBarry Smith EXTERN_C_END
4660971522cSBarry Smith 
4670971522cSBarry Smith 
4680971522cSBarry Smith EXTERN_C_BEGIN
4690971522cSBarry Smith #undef __FUNCT__
47069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
471dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
4720971522cSBarry Smith {
4730971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4740971522cSBarry Smith   PetscErrorCode    ierr;
4750971522cSBarry Smith   PetscInt          cnt = 0;
4765a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
4770971522cSBarry Smith 
4780971522cSBarry Smith   PetscFunctionBegin;
47969a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
4805a9f2f41SSatish Balay   while (ilink) {
4815a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
4825a9f2f41SSatish Balay     ilink = ilink->next;
4830971522cSBarry Smith   }
4840971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
4850971522cSBarry Smith   *n = jac->nsplits;
4860971522cSBarry Smith   PetscFunctionReturn(0);
4870971522cSBarry Smith }
4880971522cSBarry Smith EXTERN_C_END
4890971522cSBarry Smith 
4900971522cSBarry Smith #undef __FUNCT__
4910971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
4920971522cSBarry Smith /*@
4930971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
4940971522cSBarry Smith 
4950971522cSBarry Smith     Collective on PC
4960971522cSBarry Smith 
4970971522cSBarry Smith     Input Parameters:
4980971522cSBarry Smith +   pc  - the preconditioner context
4990971522cSBarry Smith .   n - the number of fields in this split
5000971522cSBarry Smith .   fields - the fields in this split
5010971522cSBarry Smith 
5020971522cSBarry Smith     Level: intermediate
5030971522cSBarry Smith 
50451f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
5050971522cSBarry Smith 
5060971522cSBarry Smith @*/
507dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
5080971522cSBarry Smith {
5090971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
5100971522cSBarry Smith 
5110971522cSBarry Smith   PetscFunctionBegin;
5120971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5130971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
5140971522cSBarry Smith   if (f) {
5150971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
5160971522cSBarry Smith   }
5170971522cSBarry Smith   PetscFunctionReturn(0);
5180971522cSBarry Smith }
5190971522cSBarry Smith 
5200971522cSBarry Smith #undef __FUNCT__
52151f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
52251f519a2SBarry Smith /*@
52351f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
52451f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
52551f519a2SBarry Smith 
52651f519a2SBarry Smith     Collective on PC
52751f519a2SBarry Smith 
52851f519a2SBarry Smith     Input Parameters:
52951f519a2SBarry Smith +   pc  - the preconditioner context
53051f519a2SBarry Smith -   bs - the block size
53151f519a2SBarry Smith 
53251f519a2SBarry Smith     Level: intermediate
53351f519a2SBarry Smith 
53451f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
53551f519a2SBarry Smith 
53651f519a2SBarry Smith @*/
53751f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
53851f519a2SBarry Smith {
53951f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
54051f519a2SBarry Smith 
54151f519a2SBarry Smith   PetscFunctionBegin;
54251f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
54351f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
54451f519a2SBarry Smith   if (f) {
54551f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
54651f519a2SBarry Smith   }
54751f519a2SBarry Smith   PetscFunctionReturn(0);
54851f519a2SBarry Smith }
54951f519a2SBarry Smith 
55051f519a2SBarry Smith #undef __FUNCT__
55169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
5520971522cSBarry Smith /*@C
55369a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
5540971522cSBarry Smith 
55569a612a9SBarry Smith    Collective on KSP
5560971522cSBarry Smith 
5570971522cSBarry Smith    Input Parameter:
5580971522cSBarry Smith .  pc - the preconditioner context
5590971522cSBarry Smith 
5600971522cSBarry Smith    Output Parameters:
5610971522cSBarry Smith +  n - the number of split
56269a612a9SBarry Smith -  pc - the array of KSP contexts
5630971522cSBarry Smith 
5640971522cSBarry Smith    Note:
56569a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
5660971522cSBarry Smith 
56769a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
5680971522cSBarry Smith 
5690971522cSBarry Smith    Level: advanced
5700971522cSBarry Smith 
5710971522cSBarry Smith .seealso: PCFIELDSPLIT
5720971522cSBarry Smith @*/
573dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
5740971522cSBarry Smith {
57569a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
5760971522cSBarry Smith 
5770971522cSBarry Smith   PetscFunctionBegin;
5780971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5790971522cSBarry Smith   PetscValidIntPointer(n,2);
58069a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
5810971522cSBarry Smith   if (f) {
58269a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
5830971522cSBarry Smith   } else {
58469a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
5850971522cSBarry Smith   }
5860971522cSBarry Smith   PetscFunctionReturn(0);
5870971522cSBarry Smith }
5880971522cSBarry Smith 
58979416396SBarry Smith EXTERN_C_BEGIN
59079416396SBarry Smith #undef __FUNCT__
59179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
592dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
59379416396SBarry Smith {
59479416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
59579416396SBarry Smith 
59679416396SBarry Smith   PetscFunctionBegin;
59779416396SBarry Smith   jac->type = type;
59879416396SBarry Smith   PetscFunctionReturn(0);
59979416396SBarry Smith }
60079416396SBarry Smith EXTERN_C_END
60179416396SBarry Smith 
60251f519a2SBarry Smith EXTERN_C_BEGIN
60351f519a2SBarry Smith #undef __FUNCT__
60451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
60551f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
60651f519a2SBarry Smith {
60751f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
60851f519a2SBarry Smith 
60951f519a2SBarry Smith   PetscFunctionBegin;
61051f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
61151f519a2SBarry 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);
61251f519a2SBarry Smith   jac->bs = bs;
61351f519a2SBarry Smith   PetscFunctionReturn(0);
61451f519a2SBarry Smith }
61551f519a2SBarry Smith EXTERN_C_END
61651f519a2SBarry Smith 
61779416396SBarry Smith #undef __FUNCT__
61879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
61979416396SBarry Smith /*@C
62079416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
62179416396SBarry Smith 
62279416396SBarry Smith    Collective on PC
62379416396SBarry Smith 
62479416396SBarry Smith    Input Parameter:
62579416396SBarry Smith .  pc - the preconditioner context
62679416396SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
62779416396SBarry Smith 
62879416396SBarry Smith    Options Database Key:
62979416396SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
63079416396SBarry Smith 
63179416396SBarry Smith    Level: Developer
63279416396SBarry Smith 
63379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
63479416396SBarry Smith 
63579416396SBarry Smith .seealso: PCCompositeSetType()
63679416396SBarry Smith 
63779416396SBarry Smith @*/
638dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
63979416396SBarry Smith {
64079416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
64179416396SBarry Smith 
64279416396SBarry Smith   PetscFunctionBegin;
64379416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
64479416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
64579416396SBarry Smith   if (f) {
64679416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
64779416396SBarry Smith   }
64879416396SBarry Smith   PetscFunctionReturn(0);
64979416396SBarry Smith }
65079416396SBarry Smith 
6510971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
6520971522cSBarry Smith /*MC
653a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
6540971522cSBarry Smith                   fields or groups of fields
6550971522cSBarry Smith 
6560971522cSBarry Smith 
6570971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
658d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
6590971522cSBarry Smith 
660a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
66169a612a9SBarry Smith          and set the options directly on the resulting KSP object
6620971522cSBarry Smith 
6630971522cSBarry Smith    Level: intermediate
6640971522cSBarry Smith 
66579416396SBarry Smith    Options Database Keys:
6662e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
6672e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
6682e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
66951f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
6702e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
67179416396SBarry Smith 
6720971522cSBarry Smith    Concepts: physics based preconditioners
6730971522cSBarry Smith 
6740971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
675c8d321e0SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
6760971522cSBarry Smith M*/
6770971522cSBarry Smith 
6780971522cSBarry Smith EXTERN_C_BEGIN
6790971522cSBarry Smith #undef __FUNCT__
6800971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
681dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
6820971522cSBarry Smith {
6830971522cSBarry Smith   PetscErrorCode ierr;
6840971522cSBarry Smith   PC_FieldSplit  *jac;
6850971522cSBarry Smith 
6860971522cSBarry Smith   PetscFunctionBegin;
6870971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
68852e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr);
6890971522cSBarry Smith   jac->bs        = -1;
6900971522cSBarry Smith   jac->nsplits   = 0;
69179416396SBarry Smith   jac->type      = PC_COMPOSITE_ADDITIVE;
69251f519a2SBarry Smith 
6930971522cSBarry Smith   pc->data     = (void*)jac;
6940971522cSBarry Smith 
6950971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
696421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
6970971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
6980971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
6990971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
7000971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
7010971522cSBarry Smith   pc->ops->applyrichardson   = 0;
7020971522cSBarry Smith 
70369a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
70469a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
7050971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
7060971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
70779416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
70879416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
70951f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
71051f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
7110971522cSBarry Smith   PetscFunctionReturn(0);
7120971522cSBarry Smith }
7130971522cSBarry Smith EXTERN_C_END
7140971522cSBarry Smith 
7150971522cSBarry Smith 
716