xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision d32f9abdbc052d6e1fd06679b17a55415c3aae30)
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;
15704ba839SBarry Smith   IS                is,cis;
16704ba839SBarry 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;
24704ba839SBarry Smith   PetscInt          nsplits;
2579416396SBarry Smith   Vec               *x,*y,w1,w2;
2697bbdb24SBarry Smith   Mat               *pmat;
27704ba839SBarry Smith   PetscTruth        issetup;
2897bbdb24SBarry Smith   PC_FieldSplitLink head;
290971522cSBarry Smith } PC_FieldSplit;
300971522cSBarry Smith 
3116913363SBarry Smith /*
3216913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
3316913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
3416913363SBarry Smith    PC you could change this.
3516913363SBarry Smith */
360971522cSBarry Smith #undef __FUNCT__
370971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
380971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
390971522cSBarry Smith {
400971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
410971522cSBarry Smith   PetscErrorCode    ierr;
420971522cSBarry Smith   PetscTruth        iascii;
430971522cSBarry Smith   PetscInt          i,j;
445a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
450971522cSBarry Smith 
460971522cSBarry Smith   PetscFunctionBegin;
470971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
480971522cSBarry Smith   if (iascii) {
4951f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
5069a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
510971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
520971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
531ab39975SBarry Smith       if (ilink->fields) {
540971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
5579416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
565a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
5779416396SBarry Smith 	  if (j > 0) {
5879416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
5979416396SBarry Smith 	  }
605a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
610971522cSBarry Smith 	}
620971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6379416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
641ab39975SBarry Smith       } else {
651ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
661ab39975SBarry Smith       }
675a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
685a9f2f41SSatish Balay       ilink = ilink->next;
690971522cSBarry Smith     }
700971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
710971522cSBarry Smith   } else {
720971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
730971522cSBarry Smith   }
740971522cSBarry Smith   PetscFunctionReturn(0);
750971522cSBarry Smith }
760971522cSBarry Smith 
770971522cSBarry Smith #undef __FUNCT__
7869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
7969a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
800971522cSBarry Smith {
810971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
820971522cSBarry Smith   PetscErrorCode    ierr;
835a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
84*d32f9abdSBarry Smith   PetscInt          i = 0,*ifields,nfields;
85*d32f9abdSBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
86*d32f9abdSBarry Smith   char              optionname[128];
870971522cSBarry Smith 
880971522cSBarry Smith   PetscFunctionBegin;
89*d32f9abdSBarry Smith   if (!ilink) {
90*d32f9abdSBarry Smith 
91521d7252SBarry Smith     if (jac->bs <= 0) {
92704ba839SBarry Smith       if (pc->pmat) {
93521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
94704ba839SBarry Smith       } else {
95704ba839SBarry Smith         jac->bs = 1;
96704ba839SBarry Smith       }
97521d7252SBarry Smith     }
98*d32f9abdSBarry Smith 
99*d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
100*d32f9abdSBarry Smith     if (!flg) {
101*d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
102*d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
103*d32f9abdSBarry Smith       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
104*d32f9abdSBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
105*d32f9abdSBarry Smith       while (PETSC_TRUE) {
106*d32f9abdSBarry Smith         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
107*d32f9abdSBarry Smith         nfields = jac->bs;
108*d32f9abdSBarry Smith         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
109*d32f9abdSBarry Smith         if (!flg2) break;
110*d32f9abdSBarry Smith         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
111*d32f9abdSBarry Smith         flg = PETSC_FALSE;
112*d32f9abdSBarry Smith         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
113*d32f9abdSBarry Smith       }
114*d32f9abdSBarry Smith       ierr = PetscFree(ifields);CHKERRQ(ierr);
115*d32f9abdSBarry Smith     }
116*d32f9abdSBarry Smith 
117*d32f9abdSBarry Smith     if (flg) {
118*d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
11979416396SBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
12079416396SBarry Smith       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
1215a9f2f41SSatish Balay       while (ilink) {
1225a9f2f41SSatish Balay 	for (i=0; i<ilink->nfields; i++) {
1235a9f2f41SSatish Balay 	  fields[ilink->fields[i]] = PETSC_TRUE;
124521d7252SBarry Smith 	}
1255a9f2f41SSatish Balay 	ilink = ilink->next;
12679416396SBarry Smith       }
12797bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
12879416396SBarry Smith       for (i=0; i<jac->bs; i++) {
12979416396SBarry Smith 	if (!fields[i]) {
13079416396SBarry Smith 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
13179416396SBarry Smith 	} else {
13279416396SBarry Smith 	  jac->defaultsplit = PETSC_FALSE;
13379416396SBarry Smith 	}
13479416396SBarry Smith       }
135521d7252SBarry Smith     }
136*d32f9abdSBarry Smith   }
13769a612a9SBarry Smith   PetscFunctionReturn(0);
13869a612a9SBarry Smith }
13969a612a9SBarry Smith 
14069a612a9SBarry Smith 
14169a612a9SBarry Smith #undef __FUNCT__
14269a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
14369a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
14469a612a9SBarry Smith {
14569a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
14669a612a9SBarry Smith   PetscErrorCode    ierr;
1475a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
148cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
14969a612a9SBarry Smith   MatStructure      flag = pc->flag;
150ccb205f8SBarry Smith   PetscTruth        sorted;
15169a612a9SBarry Smith 
15269a612a9SBarry Smith   PetscFunctionBegin;
15369a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
15497bbdb24SBarry Smith   nsplit = jac->nsplits;
1555a9f2f41SSatish Balay   ilink  = jac->head;
15697bbdb24SBarry Smith 
15797bbdb24SBarry Smith   /* get the matrices for each split */
158704ba839SBarry Smith   if (!jac->issetup) {
1591b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
16097bbdb24SBarry Smith 
161704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
162704ba839SBarry Smith 
163704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
16451f519a2SBarry Smith     bs     = jac->bs;
16597bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
166cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
1671b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1681b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1691b9fc7fcSBarry Smith       if (jac->defaultsplit) {
170704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
171704ba839SBarry Smith         ilink->csize = ccsize/nsplit;
172704ba839SBarry Smith       } else if (!ilink->is) {
173ccb205f8SBarry Smith         if (ilink->nfields > 1) {
1745a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1755a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1761b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
1771b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
1781b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
17997bbdb24SBarry Smith             }
18097bbdb24SBarry Smith           }
181704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
182ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
183ccb205f8SBarry Smith         } else {
184704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
185ccb205f8SBarry Smith         }
186704ba839SBarry Smith         ilink->csize = (ccsize/bs)*ilink->nfields;
187704ba839SBarry Smith         ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
1881b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1891b9fc7fcSBarry Smith       }
190704ba839SBarry Smith       ierr  = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
191704ba839SBarry Smith       ilink = ilink->next;
1921b9fc7fcSBarry Smith     }
1931b9fc7fcSBarry Smith   }
1941b9fc7fcSBarry Smith 
195704ba839SBarry Smith   ilink  = jac->head;
19697bbdb24SBarry Smith   if (!jac->pmat) {
197cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
198cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
199704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
200704ba839SBarry Smith       ilink = ilink->next;
201cf502942SBarry Smith     }
20297bbdb24SBarry Smith   } else {
203cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
204704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
205704ba839SBarry Smith       ilink = ilink->next;
206cf502942SBarry Smith     }
20797bbdb24SBarry Smith   }
20897bbdb24SBarry Smith 
20997bbdb24SBarry Smith   /* set up the individual PCs */
21097bbdb24SBarry Smith   i    = 0;
2115a9f2f41SSatish Balay   ilink = jac->head;
2125a9f2f41SSatish Balay   while (ilink) {
2135a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
2145a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
2155a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
21697bbdb24SBarry Smith     i++;
2175a9f2f41SSatish Balay     ilink = ilink->next;
2180971522cSBarry Smith   }
21997bbdb24SBarry Smith 
22097bbdb24SBarry Smith   /* create work vectors for each split */
2211b9fc7fcSBarry Smith   if (!jac->x) {
22279416396SBarry Smith     Vec xtmp;
22397bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
2245a9f2f41SSatish Balay     ilink = jac->head;
22597bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
226906ed7ccSBarry Smith       Vec *vl,*vr;
227906ed7ccSBarry Smith 
228906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
229906ed7ccSBarry Smith       ilink->x  = *vr;
230906ed7ccSBarry Smith       ilink->y  = *vl;
231906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
232906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
2335a9f2f41SSatish Balay       jac->x[i] = ilink->x;
2345a9f2f41SSatish Balay       jac->y[i] = ilink->y;
2355a9f2f41SSatish Balay       ilink     = ilink->next;
23697bbdb24SBarry Smith     }
23779416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
2381b9fc7fcSBarry Smith 
2395a9f2f41SSatish Balay     ilink = jac->head;
2401b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
2411b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
242704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
2435a9f2f41SSatish Balay       ilink = ilink->next;
2441b9fc7fcSBarry Smith     }
2451b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2461b9fc7fcSBarry Smith   }
2470971522cSBarry Smith   PetscFunctionReturn(0);
2480971522cSBarry Smith }
2490971522cSBarry Smith 
2505a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
251ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
252ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
2535a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
254ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
255ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
25679416396SBarry Smith 
2570971522cSBarry Smith #undef __FUNCT__
2580971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2590971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2600971522cSBarry Smith {
2610971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2620971522cSBarry Smith   PetscErrorCode    ierr;
2635a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
264e25d487eSBarry Smith   PetscInt          bs;
2650971522cSBarry Smith 
2660971522cSBarry Smith   PetscFunctionBegin;
26751f519a2SBarry Smith   CHKMEMQ;
268e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
26951f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
27051f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
27151f519a2SBarry Smith 
27279416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2731b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2740971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2755a9f2f41SSatish Balay       while (ilink) {
2765a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2775a9f2f41SSatish Balay         ilink = ilink->next;
2780971522cSBarry Smith       }
2790971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
2801b9fc7fcSBarry Smith     } else {
281efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
2825a9f2f41SSatish Balay       while (ilink) {
2835a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2845a9f2f41SSatish Balay         ilink = ilink->next;
2851b9fc7fcSBarry Smith       }
2861b9fc7fcSBarry Smith     }
28716913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
28879416396SBarry Smith     if (!jac->w1) {
28979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
29079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
29179416396SBarry Smith     }
292efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
2935a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
2945a9f2f41SSatish Balay     while (ilink->next) {
2955a9f2f41SSatish Balay       ilink = ilink->next;
2961ab39975SBarry Smith       ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
297efb30889SBarry Smith       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
2985a9f2f41SSatish Balay       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
29979416396SBarry Smith     }
30051f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
30151f519a2SBarry Smith       while (ilink->previous) {
30251f519a2SBarry Smith         ilink = ilink->previous;
3031ab39975SBarry Smith         ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
30451f519a2SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
30551f519a2SBarry Smith         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
30679416396SBarry Smith       }
30751f519a2SBarry Smith     }
30816913363SBarry Smith   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
30951f519a2SBarry Smith   CHKMEMQ;
3100971522cSBarry Smith   PetscFunctionReturn(0);
3110971522cSBarry Smith }
3120971522cSBarry Smith 
313421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
314ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
315ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
316421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
317ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
318ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
319421e10b8SBarry Smith 
320421e10b8SBarry Smith #undef __FUNCT__
321421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
322421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
323421e10b8SBarry Smith {
324421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
325421e10b8SBarry Smith   PetscErrorCode    ierr;
326421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
327421e10b8SBarry Smith   PetscInt          bs;
328421e10b8SBarry Smith 
329421e10b8SBarry Smith   PetscFunctionBegin;
330421e10b8SBarry Smith   CHKMEMQ;
331421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
332421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
333421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
334421e10b8SBarry Smith 
335421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
336421e10b8SBarry Smith     if (jac->defaultsplit) {
337421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
338421e10b8SBarry Smith       while (ilink) {
339421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
340421e10b8SBarry Smith 	ilink = ilink->next;
341421e10b8SBarry Smith       }
342421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
343421e10b8SBarry Smith     } else {
344421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
345421e10b8SBarry Smith       while (ilink) {
346421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
347421e10b8SBarry Smith 	ilink = ilink->next;
348421e10b8SBarry Smith       }
349421e10b8SBarry Smith     }
350421e10b8SBarry Smith   } else {
351421e10b8SBarry Smith     if (!jac->w1) {
352421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
353421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
354421e10b8SBarry Smith     }
355421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
356421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
357421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
358421e10b8SBarry Smith       while (ilink->next) {
359421e10b8SBarry Smith         ilink = ilink->next;
3609989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
361421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
362421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
363421e10b8SBarry Smith       }
364421e10b8SBarry Smith       while (ilink->previous) {
365421e10b8SBarry Smith         ilink = ilink->previous;
3669989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
367421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
368421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
369421e10b8SBarry Smith       }
370421e10b8SBarry Smith     } else {
371421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
372421e10b8SBarry Smith 	ilink = ilink->next;
373421e10b8SBarry Smith       }
374421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
375421e10b8SBarry Smith       while (ilink->previous) {
376421e10b8SBarry Smith 	ilink = ilink->previous;
3779989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
378421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
379421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
380421e10b8SBarry Smith       }
381421e10b8SBarry Smith     }
382421e10b8SBarry Smith   }
383421e10b8SBarry Smith   CHKMEMQ;
384421e10b8SBarry Smith   PetscFunctionReturn(0);
385421e10b8SBarry Smith }
386421e10b8SBarry Smith 
3870971522cSBarry Smith #undef __FUNCT__
3880971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
3890971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
3900971522cSBarry Smith {
3910971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3920971522cSBarry Smith   PetscErrorCode    ierr;
3935a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
3940971522cSBarry Smith 
3950971522cSBarry Smith   PetscFunctionBegin;
3965a9f2f41SSatish Balay   while (ilink) {
3975a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
3985a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
3995a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
4005a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
401704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
402704ba839SBarry Smith     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
4035a9f2f41SSatish Balay     next = ilink->next;
404704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
405704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
4065a9f2f41SSatish Balay     ilink = next;
4070971522cSBarry Smith   }
40805b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
40997bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
41079416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
41179416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
4120971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
4130971522cSBarry Smith   PetscFunctionReturn(0);
4140971522cSBarry Smith }
4150971522cSBarry Smith 
4160971522cSBarry Smith #undef __FUNCT__
4170971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
4180971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
4190971522cSBarry Smith {
4201b9fc7fcSBarry Smith   PetscErrorCode ierr;
42151f519a2SBarry Smith   PetscInt       i = 0,nfields,*fields,bs;
4221b9fc7fcSBarry Smith   PetscTruth     flg;
4231b9fc7fcSBarry Smith   char           optionname[128];
4249dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
4251b9fc7fcSBarry Smith 
4260971522cSBarry Smith   PetscFunctionBegin;
4271b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
42851f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
42951f519a2SBarry Smith   if (flg) {
43051f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
43151f519a2SBarry Smith   }
432704ba839SBarry Smith 
4339dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
434*d32f9abdSBarry Smith 
435*d32f9abdSBarry Smith   if (jac->bs > 0) {
436*d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
437*d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
43851f519a2SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
4391b9fc7fcSBarry Smith     while (PETSC_TRUE) {
44013f74950SBarry Smith       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
44151f519a2SBarry Smith       nfields = jac->bs;
4421b9fc7fcSBarry Smith       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
4431b9fc7fcSBarry Smith       if (!flg) break;
4441b9fc7fcSBarry Smith       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
4451b9fc7fcSBarry Smith       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
4461b9fc7fcSBarry Smith     }
44751f519a2SBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
448*d32f9abdSBarry Smith   }
4491b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4500971522cSBarry Smith   PetscFunctionReturn(0);
4510971522cSBarry Smith }
4520971522cSBarry Smith 
4530971522cSBarry Smith /*------------------------------------------------------------------------------------*/
4540971522cSBarry Smith 
4550971522cSBarry Smith EXTERN_C_BEGIN
4560971522cSBarry Smith #undef __FUNCT__
4570971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
458dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
4590971522cSBarry Smith {
46097bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4610971522cSBarry Smith   PetscErrorCode    ierr;
4625a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
46369a612a9SBarry Smith   char              prefix[128];
46451f519a2SBarry Smith   PetscInt          i;
4650971522cSBarry Smith 
4660971522cSBarry Smith   PetscFunctionBegin;
4670971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
46851f519a2SBarry Smith   for (i=0; i<n; i++) {
46951f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
47051f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
47151f519a2SBarry Smith   }
472704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
473704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
4745a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
4755a9f2f41SSatish Balay   ilink->nfields = n;
4765a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
4777adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
4785a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
47969a612a9SBarry Smith 
4807adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
4817adad957SLisandro Dalcin     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
48269a612a9SBarry Smith   } else {
48313f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
48469a612a9SBarry Smith   }
4855a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
4860971522cSBarry Smith 
4870971522cSBarry Smith   if (!next) {
4885a9f2f41SSatish Balay     jac->head       = ilink;
48951f519a2SBarry Smith     ilink->previous = PETSC_NULL;
4900971522cSBarry Smith   } else {
4910971522cSBarry Smith     while (next->next) {
4920971522cSBarry Smith       next = next->next;
4930971522cSBarry Smith     }
4945a9f2f41SSatish Balay     next->next      = ilink;
49551f519a2SBarry Smith     ilink->previous = next;
4960971522cSBarry Smith   }
4970971522cSBarry Smith   jac->nsplits++;
4980971522cSBarry Smith   PetscFunctionReturn(0);
4990971522cSBarry Smith }
5000971522cSBarry Smith EXTERN_C_END
5010971522cSBarry Smith 
5020971522cSBarry Smith 
5030971522cSBarry Smith EXTERN_C_BEGIN
5040971522cSBarry Smith #undef __FUNCT__
50569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
506dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
5070971522cSBarry Smith {
5080971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5090971522cSBarry Smith   PetscErrorCode    ierr;
5100971522cSBarry Smith   PetscInt          cnt = 0;
5115a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
5120971522cSBarry Smith 
5130971522cSBarry Smith   PetscFunctionBegin;
51469a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
5155a9f2f41SSatish Balay   while (ilink) {
5165a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
5175a9f2f41SSatish Balay     ilink = ilink->next;
5180971522cSBarry Smith   }
5190971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
5200971522cSBarry Smith   *n = jac->nsplits;
5210971522cSBarry Smith   PetscFunctionReturn(0);
5220971522cSBarry Smith }
5230971522cSBarry Smith EXTERN_C_END
5240971522cSBarry Smith 
525704ba839SBarry Smith EXTERN_C_BEGIN
526704ba839SBarry Smith #undef __FUNCT__
527704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
528704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
529704ba839SBarry Smith {
530704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
531704ba839SBarry Smith   PetscErrorCode    ierr;
532704ba839SBarry Smith     PC_FieldSplitLink ilink, next = jac->head;
533704ba839SBarry Smith   char              prefix[128];
534704ba839SBarry Smith 
535704ba839SBarry Smith   PetscFunctionBegin;
53616913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
5371ab39975SBarry Smith   ilink->is      = is;
5381ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
539704ba839SBarry Smith   ilink->next    = PETSC_NULL;
540704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
541704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
542704ba839SBarry Smith 
543704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
544704ba839SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
545704ba839SBarry Smith   } else {
546704ba839SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
547704ba839SBarry Smith   }
548704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
549704ba839SBarry Smith 
550704ba839SBarry Smith   if (!next) {
551704ba839SBarry Smith     jac->head       = ilink;
552704ba839SBarry Smith     ilink->previous = PETSC_NULL;
553704ba839SBarry Smith   } else {
554704ba839SBarry Smith     while (next->next) {
555704ba839SBarry Smith       next = next->next;
556704ba839SBarry Smith     }
557704ba839SBarry Smith     next->next      = ilink;
558704ba839SBarry Smith     ilink->previous = next;
559704ba839SBarry Smith   }
560704ba839SBarry Smith   jac->nsplits++;
561704ba839SBarry Smith 
562704ba839SBarry Smith   PetscFunctionReturn(0);
563704ba839SBarry Smith }
564704ba839SBarry Smith EXTERN_C_END
565704ba839SBarry Smith 
5660971522cSBarry Smith #undef __FUNCT__
5670971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
5680971522cSBarry Smith /*@
5690971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
5700971522cSBarry Smith 
5710971522cSBarry Smith     Collective on PC
5720971522cSBarry Smith 
5730971522cSBarry Smith     Input Parameters:
5740971522cSBarry Smith +   pc  - the preconditioner context
5750971522cSBarry Smith .   n - the number of fields in this split
5760971522cSBarry Smith .   fields - the fields in this split
5770971522cSBarry Smith 
5780971522cSBarry Smith     Level: intermediate
5790971522cSBarry Smith 
580*d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
581*d32f9abdSBarry Smith 
582*d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
583*d32f9abdSBarry Smith      size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean
584*d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
585*d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
586*d32f9abdSBarry Smith 
587*d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
5880971522cSBarry Smith 
5890971522cSBarry Smith @*/
590dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
5910971522cSBarry Smith {
5920971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
5930971522cSBarry Smith 
5940971522cSBarry Smith   PetscFunctionBegin;
5950971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5960971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
5970971522cSBarry Smith   if (f) {
5980971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
5990971522cSBarry Smith   }
6000971522cSBarry Smith   PetscFunctionReturn(0);
6010971522cSBarry Smith }
6020971522cSBarry Smith 
6030971522cSBarry Smith #undef __FUNCT__
604704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
605704ba839SBarry Smith /*@
606704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
607704ba839SBarry Smith 
608704ba839SBarry Smith     Collective on PC
609704ba839SBarry Smith 
610704ba839SBarry Smith     Input Parameters:
611704ba839SBarry Smith +   pc  - the preconditioner context
612704ba839SBarry Smith .   is - the index set that defines the vector elements in this field
613704ba839SBarry Smith 
614*d32f9abdSBarry Smith 
615*d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
616*d32f9abdSBarry Smith 
617704ba839SBarry Smith     Level: intermediate
618704ba839SBarry Smith 
619704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
620704ba839SBarry Smith 
621704ba839SBarry Smith @*/
622704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
623704ba839SBarry Smith {
624704ba839SBarry Smith   PetscErrorCode ierr,(*f)(PC,IS);
625704ba839SBarry Smith 
626704ba839SBarry Smith   PetscFunctionBegin;
627704ba839SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
628704ba839SBarry Smith   PetscValidHeaderSpecific(is,IS_COOKIE,1);
629704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
630704ba839SBarry Smith   if (f) {
631704ba839SBarry Smith     ierr = (*f)(pc,is);CHKERRQ(ierr);
632704ba839SBarry Smith   }
633704ba839SBarry Smith   PetscFunctionReturn(0);
634704ba839SBarry Smith }
635704ba839SBarry Smith 
636704ba839SBarry Smith #undef __FUNCT__
63751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
63851f519a2SBarry Smith /*@
63951f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
64051f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
64151f519a2SBarry Smith 
64251f519a2SBarry Smith     Collective on PC
64351f519a2SBarry Smith 
64451f519a2SBarry Smith     Input Parameters:
64551f519a2SBarry Smith +   pc  - the preconditioner context
64651f519a2SBarry Smith -   bs - the block size
64751f519a2SBarry Smith 
64851f519a2SBarry Smith     Level: intermediate
64951f519a2SBarry Smith 
65051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
65151f519a2SBarry Smith 
65251f519a2SBarry Smith @*/
65351f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
65451f519a2SBarry Smith {
65551f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
65651f519a2SBarry Smith 
65751f519a2SBarry Smith   PetscFunctionBegin;
65851f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
65951f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
66051f519a2SBarry Smith   if (f) {
66151f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
66251f519a2SBarry Smith   }
66351f519a2SBarry Smith   PetscFunctionReturn(0);
66451f519a2SBarry Smith }
66551f519a2SBarry Smith 
66651f519a2SBarry Smith #undef __FUNCT__
66769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
6680971522cSBarry Smith /*@C
66969a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
6700971522cSBarry Smith 
67169a612a9SBarry Smith    Collective on KSP
6720971522cSBarry Smith 
6730971522cSBarry Smith    Input Parameter:
6740971522cSBarry Smith .  pc - the preconditioner context
6750971522cSBarry Smith 
6760971522cSBarry Smith    Output Parameters:
6770971522cSBarry Smith +  n - the number of split
67869a612a9SBarry Smith -  pc - the array of KSP contexts
6790971522cSBarry Smith 
6800971522cSBarry Smith    Note:
681*d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
682*d32f9abdSBarry Smith    (not the KSP just the array that contains them).
6830971522cSBarry Smith 
68469a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
6850971522cSBarry Smith 
6860971522cSBarry Smith    Level: advanced
6870971522cSBarry Smith 
6880971522cSBarry Smith .seealso: PCFIELDSPLIT
6890971522cSBarry Smith @*/
690dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
6910971522cSBarry Smith {
69269a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
6930971522cSBarry Smith 
6940971522cSBarry Smith   PetscFunctionBegin;
6950971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6960971522cSBarry Smith   PetscValidIntPointer(n,2);
69769a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
6980971522cSBarry Smith   if (f) {
69969a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
7000971522cSBarry Smith   } else {
70169a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
7020971522cSBarry Smith   }
7030971522cSBarry Smith   PetscFunctionReturn(0);
7040971522cSBarry Smith }
7050971522cSBarry Smith 
70679416396SBarry Smith EXTERN_C_BEGIN
70779416396SBarry Smith #undef __FUNCT__
70879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
709dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
71079416396SBarry Smith {
71179416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
71279416396SBarry Smith 
71379416396SBarry Smith   PetscFunctionBegin;
71479416396SBarry Smith   jac->type = type;
71579416396SBarry Smith   PetscFunctionReturn(0);
71679416396SBarry Smith }
71779416396SBarry Smith EXTERN_C_END
71879416396SBarry Smith 
71951f519a2SBarry Smith EXTERN_C_BEGIN
72051f519a2SBarry Smith #undef __FUNCT__
72151f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
72251f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
72351f519a2SBarry Smith {
72451f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
72551f519a2SBarry Smith 
72651f519a2SBarry Smith   PetscFunctionBegin;
72751f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
72851f519a2SBarry 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);
72951f519a2SBarry Smith   jac->bs = bs;
73051f519a2SBarry Smith   PetscFunctionReturn(0);
73151f519a2SBarry Smith }
73251f519a2SBarry Smith EXTERN_C_END
73351f519a2SBarry Smith 
73479416396SBarry Smith #undef __FUNCT__
73579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
736bc08b0f1SBarry Smith /*@
73779416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
73879416396SBarry Smith 
73979416396SBarry Smith    Collective on PC
74079416396SBarry Smith 
74179416396SBarry Smith    Input Parameter:
74279416396SBarry Smith .  pc - the preconditioner context
743ce9499c7SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
74479416396SBarry Smith 
74579416396SBarry Smith    Options Database Key:
746ce9499c7SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
74779416396SBarry Smith 
74879416396SBarry Smith    Level: Developer
74979416396SBarry Smith 
75079416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
75179416396SBarry Smith 
75279416396SBarry Smith .seealso: PCCompositeSetType()
75379416396SBarry Smith 
75479416396SBarry Smith @*/
755dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
75679416396SBarry Smith {
75779416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
75879416396SBarry Smith 
75979416396SBarry Smith   PetscFunctionBegin;
76079416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
76179416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
76279416396SBarry Smith   if (f) {
76379416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
76479416396SBarry Smith   }
76579416396SBarry Smith   PetscFunctionReturn(0);
76679416396SBarry Smith }
76779416396SBarry Smith 
7680971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
7690971522cSBarry Smith /*MC
770a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
7710971522cSBarry Smith                   fields or groups of fields
7720971522cSBarry Smith 
7730971522cSBarry Smith 
7740971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
775d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
7760971522cSBarry Smith 
777a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
77869a612a9SBarry Smith          and set the options directly on the resulting KSP object
7790971522cSBarry Smith 
7800971522cSBarry Smith    Level: intermediate
7810971522cSBarry Smith 
78279416396SBarry Smith    Options Database Keys:
7832e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
7842e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
7852e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
78651f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
7872e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
78879416396SBarry Smith 
789*d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
790*d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
791*d32f9abdSBarry Smith 
792*d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
793*d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
794*d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
795*d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
796*d32f9abdSBarry Smith 
797*d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
798*d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
799*d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
800*d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
801*d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
802*d32f9abdSBarry Smith 
8030971522cSBarry Smith    Concepts: physics based preconditioners
8040971522cSBarry Smith 
8050971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
806*d32f9abdSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS()
8070971522cSBarry Smith M*/
8080971522cSBarry Smith 
8090971522cSBarry Smith EXTERN_C_BEGIN
8100971522cSBarry Smith #undef __FUNCT__
8110971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
812dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
8130971522cSBarry Smith {
8140971522cSBarry Smith   PetscErrorCode ierr;
8150971522cSBarry Smith   PC_FieldSplit  *jac;
8160971522cSBarry Smith 
8170971522cSBarry Smith   PetscFunctionBegin;
81838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
8190971522cSBarry Smith   jac->bs        = -1;
8200971522cSBarry Smith   jac->nsplits   = 0;
82179416396SBarry Smith   jac->type      = PC_COMPOSITE_ADDITIVE;
82251f519a2SBarry Smith 
8230971522cSBarry Smith   pc->data     = (void*)jac;
8240971522cSBarry Smith 
8250971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
826421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
8270971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
8280971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
8290971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
8300971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
8310971522cSBarry Smith   pc->ops->applyrichardson   = 0;
8320971522cSBarry Smith 
83369a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
83469a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
8350971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
8360971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
837704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
838704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
83979416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
84079416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
84151f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
84251f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
8430971522cSBarry Smith   PetscFunctionReturn(0);
8440971522cSBarry Smith }
8450971522cSBarry Smith EXTERN_C_END
8460971522cSBarry Smith 
8470971522cSBarry Smith 
848