xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 3b224e63a12ff3b822ab8c4a4b7be13110279612)
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 {
2168dd23aaSBarry 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;
2768dd23aaSBarry Smith   Mat               *Afield; /* the rows of the matrix associated with each field */
28704ba839SBarry Smith   PetscTruth        issetup;
29*3b224e63SBarry Smith   Mat               B,C,schur;   /* only used when Schur complement preconditioning is used */
30*3b224e63SBarry Smith   KSP               kspschur;
3197bbdb24SBarry Smith   PC_FieldSplitLink head;
320971522cSBarry Smith } PC_FieldSplit;
330971522cSBarry Smith 
3416913363SBarry Smith /*
3516913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
3616913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
3716913363SBarry Smith    PC you could change this.
3816913363SBarry Smith */
390971522cSBarry Smith #undef __FUNCT__
400971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
410971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
420971522cSBarry Smith {
430971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
440971522cSBarry Smith   PetscErrorCode    ierr;
450971522cSBarry Smith   PetscTruth        iascii;
460971522cSBarry Smith   PetscInt          i,j;
475a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
480971522cSBarry Smith 
490971522cSBarry Smith   PetscFunctionBegin;
500971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
510971522cSBarry Smith   if (iascii) {
5251f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
5369a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
540971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
550971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
561ab39975SBarry Smith       if (ilink->fields) {
570971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
5879416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
595a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
6079416396SBarry Smith 	  if (j > 0) {
6179416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
6279416396SBarry Smith 	  }
635a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
640971522cSBarry Smith 	}
650971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6679416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
671ab39975SBarry Smith       } else {
681ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
691ab39975SBarry Smith       }
705a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
715a9f2f41SSatish Balay       ilink = ilink->next;
720971522cSBarry Smith     }
730971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
740971522cSBarry Smith   } else {
750971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
760971522cSBarry Smith   }
770971522cSBarry Smith   PetscFunctionReturn(0);
780971522cSBarry Smith }
790971522cSBarry Smith 
800971522cSBarry Smith #undef __FUNCT__
81*3b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
82*3b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
83*3b224e63SBarry Smith {
84*3b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
85*3b224e63SBarry Smith   PetscErrorCode    ierr;
86*3b224e63SBarry Smith   PetscTruth        iascii;
87*3b224e63SBarry Smith   PetscInt          i,j;
88*3b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
89*3b224e63SBarry Smith   KSP               ksp;
90*3b224e63SBarry Smith 
91*3b224e63SBarry Smith   PetscFunctionBegin;
92*3b224e63SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
93*3b224e63SBarry Smith   if (iascii) {
94*3b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr);
95*3b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
96*3b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
97*3b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
98*3b224e63SBarry Smith       if (ilink->fields) {
99*3b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
100*3b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
101*3b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
102*3b224e63SBarry Smith 	  if (j > 0) {
103*3b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
104*3b224e63SBarry Smith 	  }
105*3b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
106*3b224e63SBarry Smith 	}
107*3b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
108*3b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
109*3b224e63SBarry Smith       } else {
110*3b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
111*3b224e63SBarry Smith       }
112*3b224e63SBarry Smith       ilink = ilink->next;
113*3b224e63SBarry Smith     }
114*3b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
115*3b224e63SBarry Smith     ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
116*3b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
117*3b224e63SBarry Smith     ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
118*3b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
119*3b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
120*3b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
121*3b224e63SBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
122*3b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
123*3b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
124*3b224e63SBarry Smith   } else {
125*3b224e63SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
126*3b224e63SBarry Smith   }
127*3b224e63SBarry Smith   PetscFunctionReturn(0);
128*3b224e63SBarry Smith }
129*3b224e63SBarry Smith 
130*3b224e63SBarry Smith #undef __FUNCT__
13169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
13269a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
1330971522cSBarry Smith {
1340971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
1350971522cSBarry Smith   PetscErrorCode    ierr;
1365a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
137d32f9abdSBarry Smith   PetscInt          i = 0,*ifields,nfields;
138d32f9abdSBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
139d32f9abdSBarry Smith   char              optionname[128];
1400971522cSBarry Smith 
1410971522cSBarry Smith   PetscFunctionBegin;
142d32f9abdSBarry Smith   if (!ilink) {
143d32f9abdSBarry Smith 
144521d7252SBarry Smith     if (jac->bs <= 0) {
145704ba839SBarry Smith       if (pc->pmat) {
146521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
147704ba839SBarry Smith       } else {
148704ba839SBarry Smith         jac->bs = 1;
149704ba839SBarry Smith       }
150521d7252SBarry Smith     }
151d32f9abdSBarry Smith 
152d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
153d32f9abdSBarry Smith     if (!flg) {
154d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
155d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
156d32f9abdSBarry Smith       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
157d32f9abdSBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
158d32f9abdSBarry Smith       while (PETSC_TRUE) {
159d32f9abdSBarry Smith         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
160d32f9abdSBarry Smith         nfields = jac->bs;
161d32f9abdSBarry Smith         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
162d32f9abdSBarry Smith         if (!flg2) break;
163d32f9abdSBarry Smith         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
164d32f9abdSBarry Smith         flg = PETSC_FALSE;
165d32f9abdSBarry Smith         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
166d32f9abdSBarry Smith       }
167d32f9abdSBarry Smith       ierr = PetscFree(ifields);CHKERRQ(ierr);
168d32f9abdSBarry Smith     }
169d32f9abdSBarry Smith 
170d32f9abdSBarry Smith     if (flg) {
171d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
17279416396SBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
17379416396SBarry Smith       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
1745a9f2f41SSatish Balay       while (ilink) {
1755a9f2f41SSatish Balay 	for (i=0; i<ilink->nfields; i++) {
1765a9f2f41SSatish Balay 	  fields[ilink->fields[i]] = PETSC_TRUE;
177521d7252SBarry Smith 	}
1785a9f2f41SSatish Balay 	ilink = ilink->next;
17979416396SBarry Smith       }
18097bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
18179416396SBarry Smith       for (i=0; i<jac->bs; i++) {
18279416396SBarry Smith 	if (!fields[i]) {
18379416396SBarry Smith 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
18479416396SBarry Smith 	} else {
18579416396SBarry Smith 	  jac->defaultsplit = PETSC_FALSE;
18679416396SBarry Smith 	}
18779416396SBarry Smith       }
18868dd23aaSBarry Smith       ierr = PetscFree(fields);CHKERRQ(ierr);
189521d7252SBarry Smith     }
190d32f9abdSBarry Smith   }
19169a612a9SBarry Smith   PetscFunctionReturn(0);
19269a612a9SBarry Smith }
19369a612a9SBarry Smith 
19469a612a9SBarry Smith 
19569a612a9SBarry Smith #undef __FUNCT__
19669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
19769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
19869a612a9SBarry Smith {
19969a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
20069a612a9SBarry Smith   PetscErrorCode    ierr;
2015a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
202cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
20369a612a9SBarry Smith   MatStructure      flag = pc->flag;
20468dd23aaSBarry Smith   PetscTruth        sorted,getsub;
20569a612a9SBarry Smith 
20669a612a9SBarry Smith   PetscFunctionBegin;
20769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
20897bbdb24SBarry Smith   nsplit = jac->nsplits;
2095a9f2f41SSatish Balay   ilink  = jac->head;
21097bbdb24SBarry Smith 
21197bbdb24SBarry Smith   /* get the matrices for each split */
212704ba839SBarry Smith   if (!jac->issetup) {
2131b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
21497bbdb24SBarry Smith 
215704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
216704ba839SBarry Smith 
217704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
21851f519a2SBarry Smith     bs     = jac->bs;
21997bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
220cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
2211b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
2221b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
2231b9fc7fcSBarry Smith       if (jac->defaultsplit) {
224704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
225704ba839SBarry Smith       } else if (!ilink->is) {
226ccb205f8SBarry Smith         if (ilink->nfields > 1) {
2275a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
2285a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
2291b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
2301b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
2311b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
23297bbdb24SBarry Smith             }
23397bbdb24SBarry Smith           }
234704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
235ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
236ccb205f8SBarry Smith         } else {
237704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
238ccb205f8SBarry Smith         }
2393e197d65SBarry Smith       }
2403e197d65SBarry Smith       ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr);
241704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
2421b9fc7fcSBarry Smith       if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
243704ba839SBarry Smith       ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
244704ba839SBarry Smith       ilink = ilink->next;
2451b9fc7fcSBarry Smith     }
2461b9fc7fcSBarry Smith   }
2471b9fc7fcSBarry Smith 
248704ba839SBarry Smith   ilink  = jac->head;
24997bbdb24SBarry Smith   if (!jac->pmat) {
250cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
251cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
252704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
253704ba839SBarry Smith       ilink = ilink->next;
254cf502942SBarry Smith     }
25597bbdb24SBarry Smith   } else {
256cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
257704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
258704ba839SBarry Smith       ilink = ilink->next;
259cf502942SBarry Smith     }
26097bbdb24SBarry Smith   }
26197bbdb24SBarry Smith 
26268dd23aaSBarry Smith   /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
26368dd23aaSBarry Smith   ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr);
264*3b224e63SBarry Smith   if (getsub && jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
26568dd23aaSBarry Smith     ilink  = jac->head;
26668dd23aaSBarry Smith     if (!jac->Afield) {
26768dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
26868dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
26911755939SBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
27068dd23aaSBarry Smith 	ilink = ilink->next;
27168dd23aaSBarry Smith       }
27268dd23aaSBarry Smith     } else {
27368dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
27411755939SBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
27568dd23aaSBarry Smith 	ilink = ilink->next;
27668dd23aaSBarry Smith       }
27768dd23aaSBarry Smith     }
27868dd23aaSBarry Smith   }
27968dd23aaSBarry Smith 
280*3b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
281*3b224e63SBarry Smith     IS       ccis;
282*3b224e63SBarry Smith     PetscInt N;
283*3b224e63SBarry Smith     if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
28468dd23aaSBarry Smith 
285*3b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
286*3b224e63SBarry Smith     if (jac->schur) {
287*3b224e63SBarry Smith       ierr  = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr);
288*3b224e63SBarry Smith       ilink = jac->head;
289*3b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
290*3b224e63SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
291*3b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
292*3b224e63SBarry Smith       ilink = ilink->next;
293*3b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
294*3b224e63SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
295*3b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
296*3b224e63SBarry Smith       ierr  = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
297*3b224e63SBarry Smith       ierr  = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,pc->flag);CHKERRQ(ierr);
298*3b224e63SBarry Smith 
299*3b224e63SBarry Smith      } else {
300*3b224e63SBarry Smith 
301*3b224e63SBarry Smith       /* extract the B and C matrices */
302*3b224e63SBarry Smith       ierr  = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr);
303*3b224e63SBarry Smith       ilink = jac->head;
304*3b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
305*3b224e63SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
306*3b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
307*3b224e63SBarry Smith       ilink = ilink->next;
308*3b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
309*3b224e63SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
310*3b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
311*3b224e63SBarry Smith       ierr  = MatCreateSchurComplement(jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr);
312*3b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
313*3b224e63SBarry Smith 
314*3b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
315*3b224e63SBarry Smith       ierr  = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
316*3b224e63SBarry Smith       ierr  = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr);
317*3b224e63SBarry Smith       ierr  = KSPAppendOptionsPrefix(jac->kspschur,"schur_");CHKERRQ(ierr);
318*3b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
319*3b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
320*3b224e63SBarry Smith 
321*3b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
322*3b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
323*3b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
324*3b224e63SBarry Smith       ilink = jac->head;
325*3b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
326*3b224e63SBarry Smith       ilink = ilink->next;
327*3b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
328*3b224e63SBarry Smith     }
329*3b224e63SBarry Smith   } else {
33097bbdb24SBarry Smith     /* set up the individual PCs */
33197bbdb24SBarry Smith     i    = 0;
3325a9f2f41SSatish Balay     ilink = jac->head;
3335a9f2f41SSatish Balay     while (ilink) {
3345a9f2f41SSatish Balay       ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
335*3b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
3365a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
3375a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
33897bbdb24SBarry Smith       i++;
3395a9f2f41SSatish Balay       ilink = ilink->next;
3400971522cSBarry Smith     }
34197bbdb24SBarry Smith 
34297bbdb24SBarry Smith     /* create work vectors for each split */
3431b9fc7fcSBarry Smith     if (!jac->x) {
34479416396SBarry Smith       Vec xtmp;
34597bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
3465a9f2f41SSatish Balay       ilink = jac->head;
34797bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
348906ed7ccSBarry Smith         Vec *vl,*vr;
349906ed7ccSBarry Smith 
350906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
351906ed7ccSBarry Smith         ilink->x  = *vr;
352906ed7ccSBarry Smith         ilink->y  = *vl;
353906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
354906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
3555a9f2f41SSatish Balay         jac->x[i] = ilink->x;
3565a9f2f41SSatish Balay         jac->y[i] = ilink->y;
3575a9f2f41SSatish Balay         ilink     = ilink->next;
35897bbdb24SBarry Smith       }
359*3b224e63SBarry Smith     }
360*3b224e63SBarry Smith   }
361*3b224e63SBarry Smith 
362*3b224e63SBarry Smith 
363*3b224e63SBarry Smith   if (1) {
364*3b224e63SBarry Smith     Vec xtmp;
365*3b224e63SBarry Smith 
36679416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
3671b9fc7fcSBarry Smith 
3685a9f2f41SSatish Balay     ilink = jac->head;
3691b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
3701b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
371704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
3725a9f2f41SSatish Balay       ilink = ilink->next;
3731b9fc7fcSBarry Smith     }
3741b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
3751b9fc7fcSBarry Smith   }
3760971522cSBarry Smith   PetscFunctionReturn(0);
3770971522cSBarry Smith }
3780971522cSBarry Smith 
3795a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
380ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
381ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
3825a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
383ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
384ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
38579416396SBarry Smith 
3860971522cSBarry Smith #undef __FUNCT__
387*3b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
388*3b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
389*3b224e63SBarry Smith {
390*3b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
391*3b224e63SBarry Smith   PetscErrorCode    ierr;
392*3b224e63SBarry Smith   KSP               ksp;
393*3b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
394*3b224e63SBarry Smith 
395*3b224e63SBarry Smith   PetscFunctionBegin;
396*3b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
397*3b224e63SBarry Smith 
398*3b224e63SBarry Smith   ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
399*3b224e63SBarry Smith   ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
400*3b224e63SBarry Smith   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
401*3b224e63SBarry Smith   ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
402*3b224e63SBarry Smith   ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
403*3b224e63SBarry Smith   ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
404*3b224e63SBarry Smith   ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
405*3b224e63SBarry Smith 
406*3b224e63SBarry Smith   ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
407*3b224e63SBarry Smith   ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
408*3b224e63SBarry Smith   ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
409*3b224e63SBarry Smith 
410*3b224e63SBarry Smith   ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
411*3b224e63SBarry Smith   ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
412*3b224e63SBarry Smith   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
413*3b224e63SBarry Smith   ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
414*3b224e63SBarry Smith   ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
415*3b224e63SBarry Smith 
416*3b224e63SBarry Smith   PetscFunctionReturn(0);
417*3b224e63SBarry Smith }
418*3b224e63SBarry Smith 
419*3b224e63SBarry Smith #undef __FUNCT__
4200971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
4210971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
4220971522cSBarry Smith {
4230971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4240971522cSBarry Smith   PetscErrorCode    ierr;
4255a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
4263e197d65SBarry Smith   PetscInt          bs,cnt;
4270971522cSBarry Smith 
4280971522cSBarry Smith   PetscFunctionBegin;
42951f519a2SBarry Smith   CHKMEMQ;
430e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
43151f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
43251f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
43351f519a2SBarry Smith 
43479416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
4351b9fc7fcSBarry Smith     if (jac->defaultsplit) {
4360971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
4375a9f2f41SSatish Balay       while (ilink) {
4385a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
4395a9f2f41SSatish Balay         ilink = ilink->next;
4400971522cSBarry Smith       }
4410971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
4421b9fc7fcSBarry Smith     } else {
443efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
4445a9f2f41SSatish Balay       while (ilink) {
4455a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
4465a9f2f41SSatish Balay         ilink = ilink->next;
4471b9fc7fcSBarry Smith       }
4481b9fc7fcSBarry Smith     }
44916913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
45079416396SBarry Smith     if (!jac->w1) {
45179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
45279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
45379416396SBarry Smith     }
454efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
4555a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
4563e197d65SBarry Smith     cnt = 1;
4575a9f2f41SSatish Balay     while (ilink->next) {
4585a9f2f41SSatish Balay       ilink = ilink->next;
4593e197d65SBarry Smith       if (jac->Afield) {
4603e197d65SBarry Smith         /* compute the residual only over the part of the vector needed */
4613e197d65SBarry Smith         ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
4623e197d65SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
4633e197d65SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4643e197d65SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4653e197d65SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
4663e197d65SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4673e197d65SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4683e197d65SBarry Smith       } else {
4693e197d65SBarry Smith         /* compute the residual over the entire vector */
4701ab39975SBarry Smith 	ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
471efb30889SBarry Smith 	ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
4725a9f2f41SSatish Balay 	ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
47379416396SBarry Smith       }
4743e197d65SBarry Smith     }
47551f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
47611755939SBarry Smith       cnt -= 2;
47751f519a2SBarry Smith       while (ilink->previous) {
47851f519a2SBarry Smith         ilink = ilink->previous;
47911755939SBarry Smith         if (jac->Afield) {
48011755939SBarry Smith 	  /* compute the residual only over the part of the vector needed */
48111755939SBarry Smith 	  ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
48211755939SBarry Smith 	  ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
48311755939SBarry Smith 	  ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
48411755939SBarry Smith 	  ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
48511755939SBarry Smith 	  ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
48611755939SBarry Smith 	  ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
48711755939SBarry Smith 	  ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
48811755939SBarry Smith         } else {
4891ab39975SBarry Smith 	  ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
49051f519a2SBarry Smith 	  ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
49151f519a2SBarry Smith 	  ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
49279416396SBarry Smith         }
49351f519a2SBarry Smith       }
49411755939SBarry Smith     }
49516913363SBarry Smith   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
49651f519a2SBarry Smith   CHKMEMQ;
4970971522cSBarry Smith   PetscFunctionReturn(0);
4980971522cSBarry Smith }
4990971522cSBarry Smith 
500421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
501ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
502ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
503421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
504ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
505ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
506421e10b8SBarry Smith 
507421e10b8SBarry Smith #undef __FUNCT__
508421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
509421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
510421e10b8SBarry Smith {
511421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
512421e10b8SBarry Smith   PetscErrorCode    ierr;
513421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
514421e10b8SBarry Smith   PetscInt          bs;
515421e10b8SBarry Smith 
516421e10b8SBarry Smith   PetscFunctionBegin;
517421e10b8SBarry Smith   CHKMEMQ;
518421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
519421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
520421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
521421e10b8SBarry Smith 
522421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
523421e10b8SBarry Smith     if (jac->defaultsplit) {
524421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
525421e10b8SBarry Smith       while (ilink) {
526421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
527421e10b8SBarry Smith 	ilink = ilink->next;
528421e10b8SBarry Smith       }
529421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
530421e10b8SBarry Smith     } else {
531421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
532421e10b8SBarry Smith       while (ilink) {
533421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
534421e10b8SBarry Smith 	ilink = ilink->next;
535421e10b8SBarry Smith       }
536421e10b8SBarry Smith     }
537421e10b8SBarry Smith   } else {
538421e10b8SBarry Smith     if (!jac->w1) {
539421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
540421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
541421e10b8SBarry Smith     }
542421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
543421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
544421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
545421e10b8SBarry Smith       while (ilink->next) {
546421e10b8SBarry Smith         ilink = ilink->next;
5479989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
548421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
549421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
550421e10b8SBarry Smith       }
551421e10b8SBarry Smith       while (ilink->previous) {
552421e10b8SBarry Smith         ilink = ilink->previous;
5539989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
554421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
555421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
556421e10b8SBarry Smith       }
557421e10b8SBarry Smith     } else {
558421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
559421e10b8SBarry Smith 	ilink = ilink->next;
560421e10b8SBarry Smith       }
561421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
562421e10b8SBarry Smith       while (ilink->previous) {
563421e10b8SBarry Smith 	ilink = ilink->previous;
5649989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
565421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
566421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
567421e10b8SBarry Smith       }
568421e10b8SBarry Smith     }
569421e10b8SBarry Smith   }
570421e10b8SBarry Smith   CHKMEMQ;
571421e10b8SBarry Smith   PetscFunctionReturn(0);
572421e10b8SBarry Smith }
573421e10b8SBarry Smith 
5740971522cSBarry Smith #undef __FUNCT__
5750971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
5760971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
5770971522cSBarry Smith {
5780971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5790971522cSBarry Smith   PetscErrorCode    ierr;
5805a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
5810971522cSBarry Smith 
5820971522cSBarry Smith   PetscFunctionBegin;
5835a9f2f41SSatish Balay   while (ilink) {
5845a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
5855a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
5865a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
5875a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
588704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
589704ba839SBarry Smith     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
5905a9f2f41SSatish Balay     next = ilink->next;
591704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
592704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
5935a9f2f41SSatish Balay     ilink = next;
5940971522cSBarry Smith   }
59505b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
59697bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
59768dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
59879416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
59979416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
600*3b224e63SBarry Smith   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
601*3b224e63SBarry Smith   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
602*3b224e63SBarry Smith   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
6030971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
6040971522cSBarry Smith   PetscFunctionReturn(0);
6050971522cSBarry Smith }
6060971522cSBarry Smith 
6070971522cSBarry Smith #undef __FUNCT__
6080971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
6090971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
6100971522cSBarry Smith {
6111b9fc7fcSBarry Smith   PetscErrorCode  ierr;
61251f519a2SBarry Smith   PetscInt        i = 0,nfields,*fields,bs;
6131b9fc7fcSBarry Smith   PetscTruth      flg;
6141b9fc7fcSBarry Smith   char            optionname[128];
6159dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
616*3b224e63SBarry Smith   PCCompositeType ctype;
6171b9fc7fcSBarry Smith 
6180971522cSBarry Smith   PetscFunctionBegin;
6191b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
62051f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
62151f519a2SBarry Smith   if (flg) {
62251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
62351f519a2SBarry Smith   }
624704ba839SBarry Smith 
625*3b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
626*3b224e63SBarry Smith   if (flg) {
627*3b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
628*3b224e63SBarry Smith   }
629d32f9abdSBarry Smith 
630d32f9abdSBarry Smith   if (jac->bs > 0) {
631d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
632d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
63351f519a2SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
6341b9fc7fcSBarry Smith     while (PETSC_TRUE) {
63513f74950SBarry Smith       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
63651f519a2SBarry Smith       nfields = jac->bs;
6371b9fc7fcSBarry Smith       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
6381b9fc7fcSBarry Smith       if (!flg) break;
6391b9fc7fcSBarry Smith       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
6401b9fc7fcSBarry Smith       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
6411b9fc7fcSBarry Smith     }
64251f519a2SBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
643d32f9abdSBarry Smith   }
6441b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6450971522cSBarry Smith   PetscFunctionReturn(0);
6460971522cSBarry Smith }
6470971522cSBarry Smith 
6480971522cSBarry Smith /*------------------------------------------------------------------------------------*/
6490971522cSBarry Smith 
6500971522cSBarry Smith EXTERN_C_BEGIN
6510971522cSBarry Smith #undef __FUNCT__
6520971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
653dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
6540971522cSBarry Smith {
65597bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6560971522cSBarry Smith   PetscErrorCode    ierr;
6575a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
65869a612a9SBarry Smith   char              prefix[128];
65951f519a2SBarry Smith   PetscInt          i;
6600971522cSBarry Smith 
6610971522cSBarry Smith   PetscFunctionBegin;
6620971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
66351f519a2SBarry Smith   for (i=0; i<n; i++) {
66451f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
66551f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
66651f519a2SBarry Smith   }
667704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
668704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
6695a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
6705a9f2f41SSatish Balay   ilink->nfields = n;
6715a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
6727adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
6735a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
67469a612a9SBarry Smith 
6757adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
6767adad957SLisandro Dalcin     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
67769a612a9SBarry Smith   } else {
67813f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
67969a612a9SBarry Smith   }
6805a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
6810971522cSBarry Smith 
6820971522cSBarry Smith   if (!next) {
6835a9f2f41SSatish Balay     jac->head       = ilink;
68451f519a2SBarry Smith     ilink->previous = PETSC_NULL;
6850971522cSBarry Smith   } else {
6860971522cSBarry Smith     while (next->next) {
6870971522cSBarry Smith       next = next->next;
6880971522cSBarry Smith     }
6895a9f2f41SSatish Balay     next->next      = ilink;
69051f519a2SBarry Smith     ilink->previous = next;
6910971522cSBarry Smith   }
6920971522cSBarry Smith   jac->nsplits++;
6930971522cSBarry Smith   PetscFunctionReturn(0);
6940971522cSBarry Smith }
6950971522cSBarry Smith EXTERN_C_END
6960971522cSBarry Smith 
6970971522cSBarry Smith 
6980971522cSBarry Smith EXTERN_C_BEGIN
6990971522cSBarry Smith #undef __FUNCT__
70069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
701dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
7020971522cSBarry Smith {
7030971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7040971522cSBarry Smith   PetscErrorCode    ierr;
7050971522cSBarry Smith   PetscInt          cnt = 0;
7065a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
7070971522cSBarry Smith 
7080971522cSBarry Smith   PetscFunctionBegin;
70969a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
7105a9f2f41SSatish Balay   while (ilink) {
7115a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
7125a9f2f41SSatish Balay     ilink = ilink->next;
7130971522cSBarry Smith   }
7140971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
7150971522cSBarry Smith   *n = jac->nsplits;
7160971522cSBarry Smith   PetscFunctionReturn(0);
7170971522cSBarry Smith }
7180971522cSBarry Smith EXTERN_C_END
7190971522cSBarry Smith 
720704ba839SBarry Smith EXTERN_C_BEGIN
721704ba839SBarry Smith #undef __FUNCT__
722704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
723704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
724704ba839SBarry Smith {
725704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
726704ba839SBarry Smith   PetscErrorCode    ierr;
727704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
728704ba839SBarry Smith   char              prefix[128];
729704ba839SBarry Smith 
730704ba839SBarry Smith   PetscFunctionBegin;
73116913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
7321ab39975SBarry Smith   ilink->is      = is;
7331ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
734704ba839SBarry Smith   ilink->next    = PETSC_NULL;
735704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
736704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
737704ba839SBarry Smith 
738704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
739704ba839SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
740704ba839SBarry Smith   } else {
741704ba839SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
742704ba839SBarry Smith   }
743704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
744704ba839SBarry Smith 
745704ba839SBarry Smith   if (!next) {
746704ba839SBarry Smith     jac->head       = ilink;
747704ba839SBarry Smith     ilink->previous = PETSC_NULL;
748704ba839SBarry Smith   } else {
749704ba839SBarry Smith     while (next->next) {
750704ba839SBarry Smith       next = next->next;
751704ba839SBarry Smith     }
752704ba839SBarry Smith     next->next      = ilink;
753704ba839SBarry Smith     ilink->previous = next;
754704ba839SBarry Smith   }
755704ba839SBarry Smith   jac->nsplits++;
756704ba839SBarry Smith 
757704ba839SBarry Smith   PetscFunctionReturn(0);
758704ba839SBarry Smith }
759704ba839SBarry Smith EXTERN_C_END
760704ba839SBarry Smith 
7610971522cSBarry Smith #undef __FUNCT__
7620971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
7630971522cSBarry Smith /*@
7640971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
7650971522cSBarry Smith 
7660971522cSBarry Smith     Collective on PC
7670971522cSBarry Smith 
7680971522cSBarry Smith     Input Parameters:
7690971522cSBarry Smith +   pc  - the preconditioner context
7700971522cSBarry Smith .   n - the number of fields in this split
7710971522cSBarry Smith .   fields - the fields in this split
7720971522cSBarry Smith 
7730971522cSBarry Smith     Level: intermediate
7740971522cSBarry Smith 
775d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
776d32f9abdSBarry Smith 
777d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
778d32f9abdSBarry 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
779d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
780d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
781d32f9abdSBarry Smith 
782d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
7830971522cSBarry Smith 
7840971522cSBarry Smith @*/
785dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
7860971522cSBarry Smith {
7870971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
7880971522cSBarry Smith 
7890971522cSBarry Smith   PetscFunctionBegin;
7900971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7910971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
7920971522cSBarry Smith   if (f) {
7930971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
7940971522cSBarry Smith   }
7950971522cSBarry Smith   PetscFunctionReturn(0);
7960971522cSBarry Smith }
7970971522cSBarry Smith 
7980971522cSBarry Smith #undef __FUNCT__
799704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
800704ba839SBarry Smith /*@
801704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
802704ba839SBarry Smith 
803704ba839SBarry Smith     Collective on PC
804704ba839SBarry Smith 
805704ba839SBarry Smith     Input Parameters:
806704ba839SBarry Smith +   pc  - the preconditioner context
807704ba839SBarry Smith .   is - the index set that defines the vector elements in this field
808704ba839SBarry Smith 
809d32f9abdSBarry Smith 
810d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
811d32f9abdSBarry Smith 
812704ba839SBarry Smith     Level: intermediate
813704ba839SBarry Smith 
814704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
815704ba839SBarry Smith 
816704ba839SBarry Smith @*/
817704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
818704ba839SBarry Smith {
819704ba839SBarry Smith   PetscErrorCode ierr,(*f)(PC,IS);
820704ba839SBarry Smith 
821704ba839SBarry Smith   PetscFunctionBegin;
822704ba839SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
823704ba839SBarry Smith   PetscValidHeaderSpecific(is,IS_COOKIE,1);
824704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
825704ba839SBarry Smith   if (f) {
826704ba839SBarry Smith     ierr = (*f)(pc,is);CHKERRQ(ierr);
827704ba839SBarry Smith   }
828704ba839SBarry Smith   PetscFunctionReturn(0);
829704ba839SBarry Smith }
830704ba839SBarry Smith 
831704ba839SBarry Smith #undef __FUNCT__
83251f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
83351f519a2SBarry Smith /*@
83451f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
83551f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
83651f519a2SBarry Smith 
83751f519a2SBarry Smith     Collective on PC
83851f519a2SBarry Smith 
83951f519a2SBarry Smith     Input Parameters:
84051f519a2SBarry Smith +   pc  - the preconditioner context
84151f519a2SBarry Smith -   bs - the block size
84251f519a2SBarry Smith 
84351f519a2SBarry Smith     Level: intermediate
84451f519a2SBarry Smith 
84551f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
84651f519a2SBarry Smith 
84751f519a2SBarry Smith @*/
84851f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
84951f519a2SBarry Smith {
85051f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
85151f519a2SBarry Smith 
85251f519a2SBarry Smith   PetscFunctionBegin;
85351f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
85451f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
85551f519a2SBarry Smith   if (f) {
85651f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
85751f519a2SBarry Smith   }
85851f519a2SBarry Smith   PetscFunctionReturn(0);
85951f519a2SBarry Smith }
86051f519a2SBarry Smith 
86151f519a2SBarry Smith #undef __FUNCT__
86269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
8630971522cSBarry Smith /*@C
86469a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
8650971522cSBarry Smith 
86669a612a9SBarry Smith    Collective on KSP
8670971522cSBarry Smith 
8680971522cSBarry Smith    Input Parameter:
8690971522cSBarry Smith .  pc - the preconditioner context
8700971522cSBarry Smith 
8710971522cSBarry Smith    Output Parameters:
8720971522cSBarry Smith +  n - the number of split
87369a612a9SBarry Smith -  pc - the array of KSP contexts
8740971522cSBarry Smith 
8750971522cSBarry Smith    Note:
876d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
877d32f9abdSBarry Smith    (not the KSP just the array that contains them).
8780971522cSBarry Smith 
87969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
8800971522cSBarry Smith 
8810971522cSBarry Smith    Level: advanced
8820971522cSBarry Smith 
8830971522cSBarry Smith .seealso: PCFIELDSPLIT
8840971522cSBarry Smith @*/
885dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
8860971522cSBarry Smith {
88769a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
8880971522cSBarry Smith 
8890971522cSBarry Smith   PetscFunctionBegin;
8900971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8910971522cSBarry Smith   PetscValidIntPointer(n,2);
89269a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
8930971522cSBarry Smith   if (f) {
89469a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
8950971522cSBarry Smith   } else {
89669a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
8970971522cSBarry Smith   }
8980971522cSBarry Smith   PetscFunctionReturn(0);
8990971522cSBarry Smith }
9000971522cSBarry Smith 
90179416396SBarry Smith EXTERN_C_BEGIN
90279416396SBarry Smith #undef __FUNCT__
90379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
904dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
90579416396SBarry Smith {
90679416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
90779416396SBarry Smith 
90879416396SBarry Smith   PetscFunctionBegin;
90979416396SBarry Smith   jac->type = type;
910*3b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
911*3b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
912*3b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
913*3b224e63SBarry Smith   } else {
914*3b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
915*3b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
916*3b224e63SBarry Smith   }
91779416396SBarry Smith   PetscFunctionReturn(0);
91879416396SBarry Smith }
91979416396SBarry Smith EXTERN_C_END
92079416396SBarry Smith 
92151f519a2SBarry Smith EXTERN_C_BEGIN
92251f519a2SBarry Smith #undef __FUNCT__
92351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
92451f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
92551f519a2SBarry Smith {
92651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
92751f519a2SBarry Smith 
92851f519a2SBarry Smith   PetscFunctionBegin;
92951f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
93051f519a2SBarry 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);
93151f519a2SBarry Smith   jac->bs = bs;
93251f519a2SBarry Smith   PetscFunctionReturn(0);
93351f519a2SBarry Smith }
93451f519a2SBarry Smith EXTERN_C_END
93551f519a2SBarry Smith 
93679416396SBarry Smith #undef __FUNCT__
93779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
938bc08b0f1SBarry Smith /*@
93979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
94079416396SBarry Smith 
94179416396SBarry Smith    Collective on PC
94279416396SBarry Smith 
94379416396SBarry Smith    Input Parameter:
94479416396SBarry Smith .  pc - the preconditioner context
945ce9499c7SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
94679416396SBarry Smith 
94779416396SBarry Smith    Options Database Key:
948ce9499c7SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
94979416396SBarry Smith 
95079416396SBarry Smith    Level: Developer
95179416396SBarry Smith 
95279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
95379416396SBarry Smith 
95479416396SBarry Smith .seealso: PCCompositeSetType()
95579416396SBarry Smith 
95679416396SBarry Smith @*/
957dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
95879416396SBarry Smith {
95979416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
96079416396SBarry Smith 
96179416396SBarry Smith   PetscFunctionBegin;
96279416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
96379416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
96479416396SBarry Smith   if (f) {
96579416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
96679416396SBarry Smith   }
96779416396SBarry Smith   PetscFunctionReturn(0);
96879416396SBarry Smith }
96979416396SBarry Smith 
9700971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
9710971522cSBarry Smith /*MC
972a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
9730971522cSBarry Smith                   fields or groups of fields
9740971522cSBarry Smith 
9750971522cSBarry Smith 
9760971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
977d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
9780971522cSBarry Smith 
979a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
98069a612a9SBarry Smith          and set the options directly on the resulting KSP object
9810971522cSBarry Smith 
9820971522cSBarry Smith    Level: intermediate
9830971522cSBarry Smith 
98479416396SBarry Smith    Options Database Keys:
9852e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
9862e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
9872e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
98851f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
9892e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
99079416396SBarry Smith 
991*3b224e63SBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -schurAblock_ and -schur,
992*3b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
993*3b224e63SBarry Smith 
994*3b224e63SBarry Smith 
995d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
996d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
997d32f9abdSBarry Smith 
998d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
999d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1000d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1001d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1002d32f9abdSBarry Smith 
1003d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
1004d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
1005d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
1006d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1007d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
1008d32f9abdSBarry Smith 
10090971522cSBarry Smith    Concepts: physics based preconditioners
10100971522cSBarry Smith 
10110971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1012d32f9abdSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS()
10130971522cSBarry Smith M*/
10140971522cSBarry Smith 
10150971522cSBarry Smith EXTERN_C_BEGIN
10160971522cSBarry Smith #undef __FUNCT__
10170971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
1018dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
10190971522cSBarry Smith {
10200971522cSBarry Smith   PetscErrorCode ierr;
10210971522cSBarry Smith   PC_FieldSplit  *jac;
10220971522cSBarry Smith 
10230971522cSBarry Smith   PetscFunctionBegin;
102438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
10250971522cSBarry Smith   jac->bs        = -1;
10260971522cSBarry Smith   jac->nsplits   = 0;
10273e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
102851f519a2SBarry Smith 
10290971522cSBarry Smith   pc->data     = (void*)jac;
10300971522cSBarry Smith 
10310971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1032421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
10330971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
10340971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
10350971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
10360971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
10370971522cSBarry Smith   pc->ops->applyrichardson   = 0;
10380971522cSBarry Smith 
103969a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
104069a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
10410971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
10420971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1043704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1044704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
104579416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
104679416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
104751f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
104851f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
10490971522cSBarry Smith   PetscFunctionReturn(0);
10500971522cSBarry Smith }
10510971522cSBarry Smith EXTERN_C_END
10520971522cSBarry Smith 
10530971522cSBarry Smith 
1054