xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 68dd23aa9c1aa4a90ae4892afbe846a842c61295)
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 {
21*68dd23aaSBarry Smith   PCCompositeType   type;
2297bbdb24SBarry Smith   PetscTruth        defaultsplit;
230971522cSBarry Smith   PetscInt          bs;
24704ba839SBarry Smith   PetscInt          nsplits;
2579416396SBarry Smith   Vec               *x,*y,w1,w2;
2697bbdb24SBarry Smith   Mat               *pmat;
27*68dd23aaSBarry Smith   Mat               *Afield; /* the rows of the matrix associated with each field */
28704ba839SBarry Smith   PetscTruth        issetup;
2997bbdb24SBarry Smith   PC_FieldSplitLink head;
300971522cSBarry Smith } PC_FieldSplit;
310971522cSBarry Smith 
3216913363SBarry Smith /*
3316913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
3416913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
3516913363SBarry Smith    PC you could change this.
3616913363SBarry Smith */
370971522cSBarry Smith #undef __FUNCT__
380971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
390971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
400971522cSBarry Smith {
410971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
420971522cSBarry Smith   PetscErrorCode    ierr;
430971522cSBarry Smith   PetscTruth        iascii;
440971522cSBarry Smith   PetscInt          i,j;
455a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
460971522cSBarry Smith 
470971522cSBarry Smith   PetscFunctionBegin;
480971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
490971522cSBarry Smith   if (iascii) {
5051f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
5169a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
520971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
530971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
541ab39975SBarry Smith       if (ilink->fields) {
550971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
5679416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
575a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
5879416396SBarry Smith 	  if (j > 0) {
5979416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
6079416396SBarry Smith 	  }
615a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
620971522cSBarry Smith 	}
630971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6479416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
651ab39975SBarry Smith       } else {
661ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
671ab39975SBarry Smith       }
685a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
695a9f2f41SSatish Balay       ilink = ilink->next;
700971522cSBarry Smith     }
710971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
720971522cSBarry Smith   } else {
730971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
740971522cSBarry Smith   }
750971522cSBarry Smith   PetscFunctionReturn(0);
760971522cSBarry Smith }
770971522cSBarry Smith 
780971522cSBarry Smith #undef __FUNCT__
7969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
8069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
810971522cSBarry Smith {
820971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
830971522cSBarry Smith   PetscErrorCode    ierr;
845a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
85d32f9abdSBarry Smith   PetscInt          i = 0,*ifields,nfields;
86d32f9abdSBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
87d32f9abdSBarry Smith   char              optionname[128];
880971522cSBarry Smith 
890971522cSBarry Smith   PetscFunctionBegin;
90d32f9abdSBarry Smith   if (!ilink) {
91d32f9abdSBarry Smith 
92521d7252SBarry Smith     if (jac->bs <= 0) {
93704ba839SBarry Smith       if (pc->pmat) {
94521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
95704ba839SBarry Smith       } else {
96704ba839SBarry Smith         jac->bs = 1;
97704ba839SBarry Smith       }
98521d7252SBarry Smith     }
99d32f9abdSBarry Smith 
100d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
101d32f9abdSBarry Smith     if (!flg) {
102d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
103d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
104d32f9abdSBarry Smith       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
105d32f9abdSBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
106d32f9abdSBarry Smith       while (PETSC_TRUE) {
107d32f9abdSBarry Smith         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
108d32f9abdSBarry Smith         nfields = jac->bs;
109d32f9abdSBarry Smith         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
110d32f9abdSBarry Smith         if (!flg2) break;
111d32f9abdSBarry Smith         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
112d32f9abdSBarry Smith         flg = PETSC_FALSE;
113d32f9abdSBarry Smith         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
114d32f9abdSBarry Smith       }
115d32f9abdSBarry Smith       ierr = PetscFree(ifields);CHKERRQ(ierr);
116d32f9abdSBarry Smith     }
117d32f9abdSBarry Smith 
118d32f9abdSBarry Smith     if (flg) {
119d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
12079416396SBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
12179416396SBarry Smith       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
1225a9f2f41SSatish Balay       while (ilink) {
1235a9f2f41SSatish Balay 	for (i=0; i<ilink->nfields; i++) {
1245a9f2f41SSatish Balay 	  fields[ilink->fields[i]] = PETSC_TRUE;
125521d7252SBarry Smith 	}
1265a9f2f41SSatish Balay 	ilink = ilink->next;
12779416396SBarry Smith       }
12897bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
12979416396SBarry Smith       for (i=0; i<jac->bs; i++) {
13079416396SBarry Smith 	if (!fields[i]) {
13179416396SBarry Smith 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
13279416396SBarry Smith 	} else {
13379416396SBarry Smith 	  jac->defaultsplit = PETSC_FALSE;
13479416396SBarry Smith 	}
13579416396SBarry Smith       }
136*68dd23aaSBarry Smith       ierr = PetscFree(fields);CHKERRQ(ierr);
137521d7252SBarry Smith     }
138d32f9abdSBarry Smith   }
13969a612a9SBarry Smith   PetscFunctionReturn(0);
14069a612a9SBarry Smith }
14169a612a9SBarry Smith 
14269a612a9SBarry Smith 
14369a612a9SBarry Smith #undef __FUNCT__
14469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
14569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
14669a612a9SBarry Smith {
14769a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
14869a612a9SBarry Smith   PetscErrorCode    ierr;
1495a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
150cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
15169a612a9SBarry Smith   MatStructure      flag = pc->flag;
152*68dd23aaSBarry Smith   PetscTruth        sorted,getsub;
15369a612a9SBarry Smith 
15469a612a9SBarry Smith   PetscFunctionBegin;
15569a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
15697bbdb24SBarry Smith   nsplit = jac->nsplits;
1575a9f2f41SSatish Balay   ilink  = jac->head;
15897bbdb24SBarry Smith 
15997bbdb24SBarry Smith   /* get the matrices for each split */
160704ba839SBarry Smith   if (!jac->issetup) {
1611b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
16297bbdb24SBarry Smith 
163704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
164704ba839SBarry Smith 
165704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
16651f519a2SBarry Smith     bs     = jac->bs;
16797bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
168cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
1691b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1701b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1711b9fc7fcSBarry Smith       if (jac->defaultsplit) {
172704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
173704ba839SBarry Smith         ilink->csize = ccsize/nsplit;
174704ba839SBarry Smith       } else if (!ilink->is) {
175ccb205f8SBarry Smith         if (ilink->nfields > 1) {
1765a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1775a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1781b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
1791b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
1801b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
18197bbdb24SBarry Smith             }
18297bbdb24SBarry Smith           }
183704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
184ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
185ccb205f8SBarry Smith         } else {
186704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
187ccb205f8SBarry Smith         }
188704ba839SBarry Smith         ilink->csize = (ccsize/bs)*ilink->nfields;
189704ba839SBarry Smith         ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
1901b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1911b9fc7fcSBarry Smith       }
192704ba839SBarry Smith       ierr  = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
193704ba839SBarry Smith       ilink = ilink->next;
1941b9fc7fcSBarry Smith     }
1951b9fc7fcSBarry Smith   }
1961b9fc7fcSBarry Smith 
197704ba839SBarry Smith   ilink  = jac->head;
19897bbdb24SBarry Smith   if (!jac->pmat) {
199cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
200cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
201704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
202704ba839SBarry Smith       ilink = ilink->next;
203cf502942SBarry Smith     }
20497bbdb24SBarry Smith   } else {
205cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
206704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
207704ba839SBarry Smith       ilink = ilink->next;
208cf502942SBarry Smith     }
20997bbdb24SBarry Smith   }
21097bbdb24SBarry Smith 
211*68dd23aaSBarry Smith   /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
212*68dd23aaSBarry Smith   ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr);
213*68dd23aaSBarry Smith   if (getsub && jac->type != PC_COMPOSITE_ADDITIVE) {
214*68dd23aaSBarry Smith     ilink  = jac->head;
215*68dd23aaSBarry Smith     if (!jac->Afield) {
216*68dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
217*68dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
218*68dd23aaSBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,ilink->csize,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
219*68dd23aaSBarry Smith 	ilink = ilink->next;
220*68dd23aaSBarry Smith       }
221*68dd23aaSBarry Smith     } else {
222*68dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
223*68dd23aaSBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,ilink->csize,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
224*68dd23aaSBarry Smith 	ilink = ilink->next;
225*68dd23aaSBarry Smith       }
226*68dd23aaSBarry Smith     }
227*68dd23aaSBarry Smith   }
228*68dd23aaSBarry Smith 
229*68dd23aaSBarry Smith 
23097bbdb24SBarry Smith   /* set up the individual PCs */
23197bbdb24SBarry Smith   i    = 0;
2325a9f2f41SSatish Balay   ilink = jac->head;
2335a9f2f41SSatish Balay   while (ilink) {
2345a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
2355a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
2365a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
23797bbdb24SBarry Smith     i++;
2385a9f2f41SSatish Balay     ilink = ilink->next;
2390971522cSBarry Smith   }
24097bbdb24SBarry Smith 
24197bbdb24SBarry Smith   /* create work vectors for each split */
2421b9fc7fcSBarry Smith   if (!jac->x) {
24379416396SBarry Smith     Vec xtmp;
24497bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
2455a9f2f41SSatish Balay     ilink = jac->head;
24697bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
247906ed7ccSBarry Smith       Vec *vl,*vr;
248906ed7ccSBarry Smith 
249906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
250906ed7ccSBarry Smith       ilink->x  = *vr;
251906ed7ccSBarry Smith       ilink->y  = *vl;
252906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
253906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
2545a9f2f41SSatish Balay       jac->x[i] = ilink->x;
2555a9f2f41SSatish Balay       jac->y[i] = ilink->y;
2565a9f2f41SSatish Balay       ilink     = ilink->next;
25797bbdb24SBarry Smith     }
25879416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
2591b9fc7fcSBarry Smith 
2605a9f2f41SSatish Balay     ilink = jac->head;
2611b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
2621b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
263704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
2645a9f2f41SSatish Balay       ilink = ilink->next;
2651b9fc7fcSBarry Smith     }
2661b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2671b9fc7fcSBarry Smith   }
2680971522cSBarry Smith   PetscFunctionReturn(0);
2690971522cSBarry Smith }
2700971522cSBarry Smith 
2715a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
272ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
273ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
2745a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
275ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
276ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
27779416396SBarry Smith 
2780971522cSBarry Smith #undef __FUNCT__
2790971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2800971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2810971522cSBarry Smith {
2820971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2830971522cSBarry Smith   PetscErrorCode    ierr;
2845a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
285e25d487eSBarry Smith   PetscInt          bs;
2860971522cSBarry Smith 
2870971522cSBarry Smith   PetscFunctionBegin;
28851f519a2SBarry Smith   CHKMEMQ;
289e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
29051f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
29151f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
29251f519a2SBarry Smith 
29379416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2941b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2950971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2965a9f2f41SSatish Balay       while (ilink) {
2975a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2985a9f2f41SSatish Balay         ilink = ilink->next;
2990971522cSBarry Smith       }
3000971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
3011b9fc7fcSBarry Smith     } else {
302efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
3035a9f2f41SSatish Balay       while (ilink) {
3045a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
3055a9f2f41SSatish Balay         ilink = ilink->next;
3061b9fc7fcSBarry Smith       }
3071b9fc7fcSBarry Smith     }
30816913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
30979416396SBarry Smith     if (!jac->w1) {
31079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
31179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
31279416396SBarry Smith     }
313efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
3145a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
3155a9f2f41SSatish Balay     while (ilink->next) {
3165a9f2f41SSatish Balay       ilink = ilink->next;
3171ab39975SBarry Smith       ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
318efb30889SBarry Smith       ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
3195a9f2f41SSatish Balay       ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
32079416396SBarry Smith     }
32151f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
32251f519a2SBarry Smith       while (ilink->previous) {
32351f519a2SBarry Smith         ilink = ilink->previous;
3241ab39975SBarry Smith         ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
32551f519a2SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
32651f519a2SBarry Smith         ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
32779416396SBarry Smith       }
32851f519a2SBarry Smith     }
32916913363SBarry Smith   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
33051f519a2SBarry Smith   CHKMEMQ;
3310971522cSBarry Smith   PetscFunctionReturn(0);
3320971522cSBarry Smith }
3330971522cSBarry Smith 
334421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
335ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
336ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
337421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
338ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
339ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
340421e10b8SBarry Smith 
341421e10b8SBarry Smith #undef __FUNCT__
342421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
343421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
344421e10b8SBarry Smith {
345421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
346421e10b8SBarry Smith   PetscErrorCode    ierr;
347421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
348421e10b8SBarry Smith   PetscInt          bs;
349421e10b8SBarry Smith 
350421e10b8SBarry Smith   PetscFunctionBegin;
351421e10b8SBarry Smith   CHKMEMQ;
352421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
353421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
354421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
355421e10b8SBarry Smith 
356421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
357421e10b8SBarry Smith     if (jac->defaultsplit) {
358421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
359421e10b8SBarry Smith       while (ilink) {
360421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
361421e10b8SBarry Smith 	ilink = ilink->next;
362421e10b8SBarry Smith       }
363421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
364421e10b8SBarry Smith     } else {
365421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
366421e10b8SBarry Smith       while (ilink) {
367421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
368421e10b8SBarry Smith 	ilink = ilink->next;
369421e10b8SBarry Smith       }
370421e10b8SBarry Smith     }
371421e10b8SBarry Smith   } else {
372421e10b8SBarry Smith     if (!jac->w1) {
373421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
374421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
375421e10b8SBarry Smith     }
376421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
377421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
378421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
379421e10b8SBarry Smith       while (ilink->next) {
380421e10b8SBarry Smith         ilink = ilink->next;
3819989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
382421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
383421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
384421e10b8SBarry Smith       }
385421e10b8SBarry Smith       while (ilink->previous) {
386421e10b8SBarry Smith         ilink = ilink->previous;
3879989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
388421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
389421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
390421e10b8SBarry Smith       }
391421e10b8SBarry Smith     } else {
392421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
393421e10b8SBarry Smith 	ilink = ilink->next;
394421e10b8SBarry Smith       }
395421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
396421e10b8SBarry Smith       while (ilink->previous) {
397421e10b8SBarry Smith 	ilink = ilink->previous;
3989989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
399421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
400421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
401421e10b8SBarry Smith       }
402421e10b8SBarry Smith     }
403421e10b8SBarry Smith   }
404421e10b8SBarry Smith   CHKMEMQ;
405421e10b8SBarry Smith   PetscFunctionReturn(0);
406421e10b8SBarry Smith }
407421e10b8SBarry Smith 
4080971522cSBarry Smith #undef __FUNCT__
4090971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
4100971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
4110971522cSBarry Smith {
4120971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4130971522cSBarry Smith   PetscErrorCode    ierr;
4145a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
4150971522cSBarry Smith 
4160971522cSBarry Smith   PetscFunctionBegin;
4175a9f2f41SSatish Balay   while (ilink) {
4185a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
4195a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
4205a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
4215a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
422704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
423704ba839SBarry Smith     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
4245a9f2f41SSatish Balay     next = ilink->next;
425704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
426704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
4275a9f2f41SSatish Balay     ilink = next;
4280971522cSBarry Smith   }
42905b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
43097bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
431*68dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
43279416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
43379416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
4340971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
4350971522cSBarry Smith   PetscFunctionReturn(0);
4360971522cSBarry Smith }
4370971522cSBarry Smith 
4380971522cSBarry Smith #undef __FUNCT__
4390971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
4400971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
4410971522cSBarry Smith {
4421b9fc7fcSBarry Smith   PetscErrorCode ierr;
44351f519a2SBarry Smith   PetscInt       i = 0,nfields,*fields,bs;
4441b9fc7fcSBarry Smith   PetscTruth     flg;
4451b9fc7fcSBarry Smith   char           optionname[128];
4469dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
4471b9fc7fcSBarry Smith 
4480971522cSBarry Smith   PetscFunctionBegin;
4491b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
45051f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
45151f519a2SBarry Smith   if (flg) {
45251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
45351f519a2SBarry Smith   }
454704ba839SBarry Smith 
4559dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
456d32f9abdSBarry Smith 
457d32f9abdSBarry Smith   if (jac->bs > 0) {
458d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
459d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
46051f519a2SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
4611b9fc7fcSBarry Smith     while (PETSC_TRUE) {
46213f74950SBarry Smith       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
46351f519a2SBarry Smith       nfields = jac->bs;
4641b9fc7fcSBarry Smith       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
4651b9fc7fcSBarry Smith       if (!flg) break;
4661b9fc7fcSBarry Smith       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
4671b9fc7fcSBarry Smith       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
4681b9fc7fcSBarry Smith     }
46951f519a2SBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
470d32f9abdSBarry Smith   }
4711b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4720971522cSBarry Smith   PetscFunctionReturn(0);
4730971522cSBarry Smith }
4740971522cSBarry Smith 
4750971522cSBarry Smith /*------------------------------------------------------------------------------------*/
4760971522cSBarry Smith 
4770971522cSBarry Smith EXTERN_C_BEGIN
4780971522cSBarry Smith #undef __FUNCT__
4790971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
480dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
4810971522cSBarry Smith {
48297bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4830971522cSBarry Smith   PetscErrorCode    ierr;
4845a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
48569a612a9SBarry Smith   char              prefix[128];
48651f519a2SBarry Smith   PetscInt          i;
4870971522cSBarry Smith 
4880971522cSBarry Smith   PetscFunctionBegin;
4890971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
49051f519a2SBarry Smith   for (i=0; i<n; i++) {
49151f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
49251f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
49351f519a2SBarry Smith   }
494704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
495704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
4965a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
4975a9f2f41SSatish Balay   ilink->nfields = n;
4985a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
4997adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
5005a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
50169a612a9SBarry Smith 
5027adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
5037adad957SLisandro Dalcin     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
50469a612a9SBarry Smith   } else {
50513f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
50669a612a9SBarry Smith   }
5075a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
5080971522cSBarry Smith 
5090971522cSBarry Smith   if (!next) {
5105a9f2f41SSatish Balay     jac->head       = ilink;
51151f519a2SBarry Smith     ilink->previous = PETSC_NULL;
5120971522cSBarry Smith   } else {
5130971522cSBarry Smith     while (next->next) {
5140971522cSBarry Smith       next = next->next;
5150971522cSBarry Smith     }
5165a9f2f41SSatish Balay     next->next      = ilink;
51751f519a2SBarry Smith     ilink->previous = next;
5180971522cSBarry Smith   }
5190971522cSBarry Smith   jac->nsplits++;
5200971522cSBarry Smith   PetscFunctionReturn(0);
5210971522cSBarry Smith }
5220971522cSBarry Smith EXTERN_C_END
5230971522cSBarry Smith 
5240971522cSBarry Smith 
5250971522cSBarry Smith EXTERN_C_BEGIN
5260971522cSBarry Smith #undef __FUNCT__
52769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
528dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
5290971522cSBarry Smith {
5300971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5310971522cSBarry Smith   PetscErrorCode    ierr;
5320971522cSBarry Smith   PetscInt          cnt = 0;
5335a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
5340971522cSBarry Smith 
5350971522cSBarry Smith   PetscFunctionBegin;
53669a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
5375a9f2f41SSatish Balay   while (ilink) {
5385a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
5395a9f2f41SSatish Balay     ilink = ilink->next;
5400971522cSBarry Smith   }
5410971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
5420971522cSBarry Smith   *n = jac->nsplits;
5430971522cSBarry Smith   PetscFunctionReturn(0);
5440971522cSBarry Smith }
5450971522cSBarry Smith EXTERN_C_END
5460971522cSBarry Smith 
547704ba839SBarry Smith EXTERN_C_BEGIN
548704ba839SBarry Smith #undef __FUNCT__
549704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
550704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
551704ba839SBarry Smith {
552704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
553704ba839SBarry Smith   PetscErrorCode    ierr;
554704ba839SBarry Smith     PC_FieldSplitLink ilink, next = jac->head;
555704ba839SBarry Smith   char              prefix[128];
556704ba839SBarry Smith 
557704ba839SBarry Smith   PetscFunctionBegin;
55816913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
5591ab39975SBarry Smith   ilink->is      = is;
5601ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
561704ba839SBarry Smith   ilink->next    = PETSC_NULL;
562704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
563704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
564704ba839SBarry Smith 
565704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
566704ba839SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
567704ba839SBarry Smith   } else {
568704ba839SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
569704ba839SBarry Smith   }
570704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
571704ba839SBarry Smith 
572704ba839SBarry Smith   if (!next) {
573704ba839SBarry Smith     jac->head       = ilink;
574704ba839SBarry Smith     ilink->previous = PETSC_NULL;
575704ba839SBarry Smith   } else {
576704ba839SBarry Smith     while (next->next) {
577704ba839SBarry Smith       next = next->next;
578704ba839SBarry Smith     }
579704ba839SBarry Smith     next->next      = ilink;
580704ba839SBarry Smith     ilink->previous = next;
581704ba839SBarry Smith   }
582704ba839SBarry Smith   jac->nsplits++;
583704ba839SBarry Smith 
584704ba839SBarry Smith   PetscFunctionReturn(0);
585704ba839SBarry Smith }
586704ba839SBarry Smith EXTERN_C_END
587704ba839SBarry Smith 
5880971522cSBarry Smith #undef __FUNCT__
5890971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
5900971522cSBarry Smith /*@
5910971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
5920971522cSBarry Smith 
5930971522cSBarry Smith     Collective on PC
5940971522cSBarry Smith 
5950971522cSBarry Smith     Input Parameters:
5960971522cSBarry Smith +   pc  - the preconditioner context
5970971522cSBarry Smith .   n - the number of fields in this split
5980971522cSBarry Smith .   fields - the fields in this split
5990971522cSBarry Smith 
6000971522cSBarry Smith     Level: intermediate
6010971522cSBarry Smith 
602d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
603d32f9abdSBarry Smith 
604d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
605d32f9abdSBarry 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
606d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
607d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
608d32f9abdSBarry Smith 
609d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
6100971522cSBarry Smith 
6110971522cSBarry Smith @*/
612dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
6130971522cSBarry Smith {
6140971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
6150971522cSBarry Smith 
6160971522cSBarry Smith   PetscFunctionBegin;
6170971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6180971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
6190971522cSBarry Smith   if (f) {
6200971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
6210971522cSBarry Smith   }
6220971522cSBarry Smith   PetscFunctionReturn(0);
6230971522cSBarry Smith }
6240971522cSBarry Smith 
6250971522cSBarry Smith #undef __FUNCT__
626704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
627704ba839SBarry Smith /*@
628704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
629704ba839SBarry Smith 
630704ba839SBarry Smith     Collective on PC
631704ba839SBarry Smith 
632704ba839SBarry Smith     Input Parameters:
633704ba839SBarry Smith +   pc  - the preconditioner context
634704ba839SBarry Smith .   is - the index set that defines the vector elements in this field
635704ba839SBarry Smith 
636d32f9abdSBarry Smith 
637d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
638d32f9abdSBarry Smith 
639704ba839SBarry Smith     Level: intermediate
640704ba839SBarry Smith 
641704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
642704ba839SBarry Smith 
643704ba839SBarry Smith @*/
644704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
645704ba839SBarry Smith {
646704ba839SBarry Smith   PetscErrorCode ierr,(*f)(PC,IS);
647704ba839SBarry Smith 
648704ba839SBarry Smith   PetscFunctionBegin;
649704ba839SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
650704ba839SBarry Smith   PetscValidHeaderSpecific(is,IS_COOKIE,1);
651704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
652704ba839SBarry Smith   if (f) {
653704ba839SBarry Smith     ierr = (*f)(pc,is);CHKERRQ(ierr);
654704ba839SBarry Smith   }
655704ba839SBarry Smith   PetscFunctionReturn(0);
656704ba839SBarry Smith }
657704ba839SBarry Smith 
658704ba839SBarry Smith #undef __FUNCT__
65951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
66051f519a2SBarry Smith /*@
66151f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
66251f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
66351f519a2SBarry Smith 
66451f519a2SBarry Smith     Collective on PC
66551f519a2SBarry Smith 
66651f519a2SBarry Smith     Input Parameters:
66751f519a2SBarry Smith +   pc  - the preconditioner context
66851f519a2SBarry Smith -   bs - the block size
66951f519a2SBarry Smith 
67051f519a2SBarry Smith     Level: intermediate
67151f519a2SBarry Smith 
67251f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
67351f519a2SBarry Smith 
67451f519a2SBarry Smith @*/
67551f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
67651f519a2SBarry Smith {
67751f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
67851f519a2SBarry Smith 
67951f519a2SBarry Smith   PetscFunctionBegin;
68051f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
68151f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
68251f519a2SBarry Smith   if (f) {
68351f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
68451f519a2SBarry Smith   }
68551f519a2SBarry Smith   PetscFunctionReturn(0);
68651f519a2SBarry Smith }
68751f519a2SBarry Smith 
68851f519a2SBarry Smith #undef __FUNCT__
68969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
6900971522cSBarry Smith /*@C
69169a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
6920971522cSBarry Smith 
69369a612a9SBarry Smith    Collective on KSP
6940971522cSBarry Smith 
6950971522cSBarry Smith    Input Parameter:
6960971522cSBarry Smith .  pc - the preconditioner context
6970971522cSBarry Smith 
6980971522cSBarry Smith    Output Parameters:
6990971522cSBarry Smith +  n - the number of split
70069a612a9SBarry Smith -  pc - the array of KSP contexts
7010971522cSBarry Smith 
7020971522cSBarry Smith    Note:
703d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
704d32f9abdSBarry Smith    (not the KSP just the array that contains them).
7050971522cSBarry Smith 
70669a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
7070971522cSBarry Smith 
7080971522cSBarry Smith    Level: advanced
7090971522cSBarry Smith 
7100971522cSBarry Smith .seealso: PCFIELDSPLIT
7110971522cSBarry Smith @*/
712dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
7130971522cSBarry Smith {
71469a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
7150971522cSBarry Smith 
7160971522cSBarry Smith   PetscFunctionBegin;
7170971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7180971522cSBarry Smith   PetscValidIntPointer(n,2);
71969a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
7200971522cSBarry Smith   if (f) {
72169a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
7220971522cSBarry Smith   } else {
72369a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
7240971522cSBarry Smith   }
7250971522cSBarry Smith   PetscFunctionReturn(0);
7260971522cSBarry Smith }
7270971522cSBarry Smith 
72879416396SBarry Smith EXTERN_C_BEGIN
72979416396SBarry Smith #undef __FUNCT__
73079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
731dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
73279416396SBarry Smith {
73379416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
73479416396SBarry Smith 
73579416396SBarry Smith   PetscFunctionBegin;
73679416396SBarry Smith   jac->type = type;
73779416396SBarry Smith   PetscFunctionReturn(0);
73879416396SBarry Smith }
73979416396SBarry Smith EXTERN_C_END
74079416396SBarry Smith 
74151f519a2SBarry Smith EXTERN_C_BEGIN
74251f519a2SBarry Smith #undef __FUNCT__
74351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
74451f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
74551f519a2SBarry Smith {
74651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
74751f519a2SBarry Smith 
74851f519a2SBarry Smith   PetscFunctionBegin;
74951f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
75051f519a2SBarry 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);
75151f519a2SBarry Smith   jac->bs = bs;
75251f519a2SBarry Smith   PetscFunctionReturn(0);
75351f519a2SBarry Smith }
75451f519a2SBarry Smith EXTERN_C_END
75551f519a2SBarry Smith 
75679416396SBarry Smith #undef __FUNCT__
75779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
758bc08b0f1SBarry Smith /*@
75979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
76079416396SBarry Smith 
76179416396SBarry Smith    Collective on PC
76279416396SBarry Smith 
76379416396SBarry Smith    Input Parameter:
76479416396SBarry Smith .  pc - the preconditioner context
765ce9499c7SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
76679416396SBarry Smith 
76779416396SBarry Smith    Options Database Key:
768ce9499c7SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
76979416396SBarry Smith 
77079416396SBarry Smith    Level: Developer
77179416396SBarry Smith 
77279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
77379416396SBarry Smith 
77479416396SBarry Smith .seealso: PCCompositeSetType()
77579416396SBarry Smith 
77679416396SBarry Smith @*/
777dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
77879416396SBarry Smith {
77979416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
78079416396SBarry Smith 
78179416396SBarry Smith   PetscFunctionBegin;
78279416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
78379416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
78479416396SBarry Smith   if (f) {
78579416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
78679416396SBarry Smith   }
78779416396SBarry Smith   PetscFunctionReturn(0);
78879416396SBarry Smith }
78979416396SBarry Smith 
7900971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
7910971522cSBarry Smith /*MC
792a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
7930971522cSBarry Smith                   fields or groups of fields
7940971522cSBarry Smith 
7950971522cSBarry Smith 
7960971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
797d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
7980971522cSBarry Smith 
799a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
80069a612a9SBarry Smith          and set the options directly on the resulting KSP object
8010971522cSBarry Smith 
8020971522cSBarry Smith    Level: intermediate
8030971522cSBarry Smith 
80479416396SBarry Smith    Options Database Keys:
8052e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
8062e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
8072e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
80851f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
8092e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
81079416396SBarry Smith 
811d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
812d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
813d32f9abdSBarry Smith 
814d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
815d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
816d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
817d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
818d32f9abdSBarry Smith 
819d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
820d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
821d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
822d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
823d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
824d32f9abdSBarry Smith 
8250971522cSBarry Smith    Concepts: physics based preconditioners
8260971522cSBarry Smith 
8270971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
828d32f9abdSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS()
8290971522cSBarry Smith M*/
8300971522cSBarry Smith 
8310971522cSBarry Smith EXTERN_C_BEGIN
8320971522cSBarry Smith #undef __FUNCT__
8330971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
834dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
8350971522cSBarry Smith {
8360971522cSBarry Smith   PetscErrorCode ierr;
8370971522cSBarry Smith   PC_FieldSplit  *jac;
8380971522cSBarry Smith 
8390971522cSBarry Smith   PetscFunctionBegin;
84038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
8410971522cSBarry Smith   jac->bs        = -1;
8420971522cSBarry Smith   jac->nsplits   = 0;
84379416396SBarry Smith   jac->type      = PC_COMPOSITE_ADDITIVE;
84451f519a2SBarry Smith 
8450971522cSBarry Smith   pc->data     = (void*)jac;
8460971522cSBarry Smith 
8470971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
848421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
8490971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
8500971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
8510971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
8520971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
8530971522cSBarry Smith   pc->ops->applyrichardson   = 0;
8540971522cSBarry Smith 
85569a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
85669a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
8570971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
8580971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
859704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
860704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
86179416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
86279416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
86351f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
86451f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
8650971522cSBarry Smith   PetscFunctionReturn(0);
8660971522cSBarry Smith }
8670971522cSBarry Smith EXTERN_C_END
8680971522cSBarry Smith 
8690971522cSBarry Smith 
870