xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 9e7d6b0a155aa3071c171c8ef5a75198ca92b5e7)
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;
293b224e63SBarry Smith   Mat               B,C,schur;   /* only used when Schur complement preconditioning is used */
303b224e63SBarry Smith   KSP               kspschur;
31e69d4d44SBarry Smith   PetscTruth        schurpre;    /* preconditioner for the Schur complement is built from pmat[1] == D */
3297bbdb24SBarry Smith   PC_FieldSplitLink head;
330971522cSBarry Smith } PC_FieldSplit;
340971522cSBarry Smith 
3516913363SBarry Smith /*
3616913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
3716913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
3816913363SBarry Smith    PC you could change this.
3916913363SBarry Smith */
400971522cSBarry Smith #undef __FUNCT__
410971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
420971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
430971522cSBarry Smith {
440971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
450971522cSBarry Smith   PetscErrorCode    ierr;
460971522cSBarry Smith   PetscTruth        iascii;
470971522cSBarry Smith   PetscInt          i,j;
485a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
490971522cSBarry Smith 
500971522cSBarry Smith   PetscFunctionBegin;
510971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
520971522cSBarry Smith   if (iascii) {
5351f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
5469a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
550971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
560971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
571ab39975SBarry Smith       if (ilink->fields) {
580971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
5979416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
605a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
6179416396SBarry Smith 	  if (j > 0) {
6279416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
6379416396SBarry Smith 	  }
645a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
650971522cSBarry Smith 	}
660971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6779416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
681ab39975SBarry Smith       } else {
691ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
701ab39975SBarry Smith       }
715a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
725a9f2f41SSatish Balay       ilink = ilink->next;
730971522cSBarry Smith     }
740971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
750971522cSBarry Smith   } else {
760971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
770971522cSBarry Smith   }
780971522cSBarry Smith   PetscFunctionReturn(0);
790971522cSBarry Smith }
800971522cSBarry Smith 
810971522cSBarry Smith #undef __FUNCT__
823b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
833b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
843b224e63SBarry Smith {
853b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
863b224e63SBarry Smith   PetscErrorCode    ierr;
873b224e63SBarry Smith   PetscTruth        iascii;
883b224e63SBarry Smith   PetscInt          i,j;
893b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
903b224e63SBarry Smith   KSP               ksp;
913b224e63SBarry Smith 
923b224e63SBarry Smith   PetscFunctionBegin;
933b224e63SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
943b224e63SBarry Smith   if (iascii) {
953b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr);
963b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
973b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
983b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
993b224e63SBarry Smith       if (ilink->fields) {
1003b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1013b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1023b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1033b224e63SBarry Smith 	  if (j > 0) {
1043b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1053b224e63SBarry Smith 	  }
1063b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1073b224e63SBarry Smith 	}
1083b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1093b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1103b224e63SBarry Smith       } else {
1113b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1123b224e63SBarry Smith       }
1133b224e63SBarry Smith       ilink = ilink->next;
1143b224e63SBarry Smith     }
1153b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
1163b224e63SBarry Smith     ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1173b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1183b224e63SBarry Smith     ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
1193b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1203b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
1213b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1223b224e63SBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
1233b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1243b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1253b224e63SBarry Smith   } else {
1263b224e63SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1273b224e63SBarry Smith   }
1283b224e63SBarry Smith   PetscFunctionReturn(0);
1293b224e63SBarry Smith }
1303b224e63SBarry Smith 
1313b224e63SBarry Smith #undef __FUNCT__
13269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
13369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
1340971522cSBarry Smith {
1350971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
1360971522cSBarry Smith   PetscErrorCode    ierr;
1375a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
138d32f9abdSBarry Smith   PetscInt          i = 0,*ifields,nfields;
139d32f9abdSBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
140d32f9abdSBarry Smith   char              optionname[128];
1410971522cSBarry Smith 
1420971522cSBarry Smith   PetscFunctionBegin;
143d32f9abdSBarry Smith   if (!ilink) {
144d32f9abdSBarry Smith 
145521d7252SBarry Smith     if (jac->bs <= 0) {
146704ba839SBarry Smith       if (pc->pmat) {
147521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
148704ba839SBarry Smith       } else {
149704ba839SBarry Smith         jac->bs = 1;
150704ba839SBarry Smith       }
151521d7252SBarry Smith     }
152d32f9abdSBarry Smith 
153d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
154d32f9abdSBarry Smith     if (!flg) {
155d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
156d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
157d32f9abdSBarry Smith       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
158d32f9abdSBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
159d32f9abdSBarry Smith       while (PETSC_TRUE) {
160d32f9abdSBarry Smith         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
161d32f9abdSBarry Smith         nfields = jac->bs;
162d32f9abdSBarry Smith         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
163d32f9abdSBarry Smith         if (!flg2) break;
164d32f9abdSBarry Smith         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
165d32f9abdSBarry Smith         flg = PETSC_FALSE;
166d32f9abdSBarry Smith         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
167d32f9abdSBarry Smith       }
168d32f9abdSBarry Smith       ierr = PetscFree(ifields);CHKERRQ(ierr);
169d32f9abdSBarry Smith     }
170d32f9abdSBarry Smith 
171d32f9abdSBarry Smith     if (flg) {
172d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
17379416396SBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
17479416396SBarry Smith       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
1755a9f2f41SSatish Balay       while (ilink) {
1765a9f2f41SSatish Balay 	for (i=0; i<ilink->nfields; i++) {
1775a9f2f41SSatish Balay 	  fields[ilink->fields[i]] = PETSC_TRUE;
178521d7252SBarry Smith 	}
1795a9f2f41SSatish Balay 	ilink = ilink->next;
18079416396SBarry Smith       }
18197bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
18279416396SBarry Smith       for (i=0; i<jac->bs; i++) {
18379416396SBarry Smith 	if (!fields[i]) {
18479416396SBarry Smith 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
18579416396SBarry Smith 	} else {
18679416396SBarry Smith 	  jac->defaultsplit = PETSC_FALSE;
18779416396SBarry Smith 	}
18879416396SBarry Smith       }
18968dd23aaSBarry Smith       ierr = PetscFree(fields);CHKERRQ(ierr);
190521d7252SBarry Smith     }
191d32f9abdSBarry Smith   }
19269a612a9SBarry Smith   PetscFunctionReturn(0);
19369a612a9SBarry Smith }
19469a612a9SBarry Smith 
19569a612a9SBarry Smith 
19669a612a9SBarry Smith #undef __FUNCT__
19769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
19869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
19969a612a9SBarry Smith {
20069a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
20169a612a9SBarry Smith   PetscErrorCode    ierr;
2025a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
203cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
20469a612a9SBarry Smith   MatStructure      flag = pc->flag;
20568dd23aaSBarry Smith   PetscTruth        sorted,getsub;
20669a612a9SBarry Smith 
20769a612a9SBarry Smith   PetscFunctionBegin;
20869a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
20997bbdb24SBarry Smith   nsplit = jac->nsplits;
2105a9f2f41SSatish Balay   ilink  = jac->head;
21197bbdb24SBarry Smith 
21297bbdb24SBarry Smith   /* get the matrices for each split */
213704ba839SBarry Smith   if (!jac->issetup) {
2141b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
21597bbdb24SBarry Smith 
216704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
217704ba839SBarry Smith 
218704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
21951f519a2SBarry Smith     bs     = jac->bs;
22097bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
221cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
2221b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
2231b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
2241b9fc7fcSBarry Smith       if (jac->defaultsplit) {
225704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
226704ba839SBarry Smith       } else if (!ilink->is) {
227ccb205f8SBarry Smith         if (ilink->nfields > 1) {
2285a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
2295a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
2301b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
2311b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
2321b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
23397bbdb24SBarry Smith             }
23497bbdb24SBarry Smith           }
235704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
236ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
237ccb205f8SBarry Smith         } else {
238704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
239ccb205f8SBarry Smith         }
2403e197d65SBarry Smith       }
2413e197d65SBarry Smith       ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr);
242704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
2431b9fc7fcSBarry Smith       if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
244704ba839SBarry Smith       ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
245704ba839SBarry Smith       ilink = ilink->next;
2461b9fc7fcSBarry Smith     }
2471b9fc7fcSBarry Smith   }
2481b9fc7fcSBarry Smith 
249704ba839SBarry Smith   ilink  = jac->head;
25097bbdb24SBarry Smith   if (!jac->pmat) {
251cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
252cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
253704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
254704ba839SBarry Smith       ilink = ilink->next;
255cf502942SBarry Smith     }
25697bbdb24SBarry Smith   } else {
257cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
258704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
259704ba839SBarry Smith       ilink = ilink->next;
260cf502942SBarry Smith     }
26197bbdb24SBarry Smith   }
26297bbdb24SBarry Smith 
26368dd23aaSBarry Smith   /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
26468dd23aaSBarry Smith   ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr);
2653b224e63SBarry Smith   if (getsub && jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
26668dd23aaSBarry Smith     ilink  = jac->head;
26768dd23aaSBarry Smith     if (!jac->Afield) {
26868dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
26968dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
27011755939SBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
27168dd23aaSBarry Smith 	ilink = ilink->next;
27268dd23aaSBarry Smith       }
27368dd23aaSBarry Smith     } else {
27468dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
27511755939SBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
27668dd23aaSBarry Smith 	ilink = ilink->next;
27768dd23aaSBarry Smith       }
27868dd23aaSBarry Smith     }
27968dd23aaSBarry Smith   }
28068dd23aaSBarry Smith 
2813b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
2823b224e63SBarry Smith     IS       ccis;
283e69d4d44SBarry Smith     PetscInt N,nlocal,nis;
2843b224e63SBarry Smith     if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
28568dd23aaSBarry Smith 
2863b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
2873b224e63SBarry Smith     if (jac->schur) {
2883b224e63SBarry Smith       ierr  = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr);
2893b224e63SBarry Smith       ilink = jac->head;
2903b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
291e69d4d44SBarry Smith       ierr  = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr);
292e69d4d44SBarry Smith       ierr  = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr);
293e69d4d44SBarry Smith       nlocal = nlocal - nis;
294e69d4d44SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
2953b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
2963b224e63SBarry Smith       ilink = ilink->next;
2973b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
298e69d4d44SBarry Smith       ierr  = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr);
299e69d4d44SBarry Smith       ierr  = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr);
300e69d4d44SBarry Smith       nlocal = nlocal - nis;
301e69d4d44SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
3023b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3033b224e63SBarry Smith       ierr  = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
304e69d4d44SBarry Smith       if (jac->schurpre) {
305e69d4d44SBarry Smith         ierr  = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],pc->flag);CHKERRQ(ierr);
306e69d4d44SBarry Smith       } else {
3073b224e63SBarry Smith         ierr  = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,pc->flag);CHKERRQ(ierr);
308e69d4d44SBarry Smith       }
3093b224e63SBarry Smith 
3103b224e63SBarry Smith      } else {
3111cee3971SBarry Smith       KSP ksp;
3123b224e63SBarry Smith 
3133b224e63SBarry Smith       /* extract the B and C matrices */
3143b224e63SBarry Smith       ierr  = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr);
3153b224e63SBarry Smith       ilink = jac->head;
3163b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
317e69d4d44SBarry Smith       ierr  = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr);
318e69d4d44SBarry Smith       ierr  = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr);
319e69d4d44SBarry Smith       nlocal = nlocal - nis;
320e69d4d44SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
3213b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3223b224e63SBarry Smith       ilink = ilink->next;
3233b224e63SBarry Smith       ierr  = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr);
324e69d4d44SBarry Smith       ierr  = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr);
325e69d4d44SBarry Smith       ierr  = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr);
326e69d4d44SBarry Smith       nlocal = nlocal - nis;
327e69d4d44SBarry Smith       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
3283b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3293b224e63SBarry Smith       ierr  = MatCreateSchurComplement(jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr);
3301cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
331e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
3323b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
3333b224e63SBarry Smith 
3343b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
3351cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
336e69d4d44SBarry Smith       if (jac->schurpre) {
337e69d4d44SBarry Smith         ierr  = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
338e69d4d44SBarry Smith       } else {
3393b224e63SBarry Smith         ierr  = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
340e69d4d44SBarry Smith       }
3413b224e63SBarry Smith       ierr  = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr);
3423b224e63SBarry Smith       ierr  = KSPAppendOptionsPrefix(jac->kspschur,"schur_");CHKERRQ(ierr);
3433b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
3443b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
3453b224e63SBarry Smith 
3463b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
3473b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
3483b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
3493b224e63SBarry Smith       ilink = jac->head;
3503b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
3513b224e63SBarry Smith       ilink = ilink->next;
3523b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
3533b224e63SBarry Smith     }
3543b224e63SBarry Smith   } else {
35597bbdb24SBarry Smith     /* set up the individual PCs */
35697bbdb24SBarry Smith     i    = 0;
3575a9f2f41SSatish Balay     ilink = jac->head;
3585a9f2f41SSatish Balay     while (ilink) {
3595a9f2f41SSatish Balay       ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
3603b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
3615a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
3625a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
36397bbdb24SBarry Smith       i++;
3645a9f2f41SSatish Balay       ilink = ilink->next;
3650971522cSBarry Smith     }
36697bbdb24SBarry Smith 
36797bbdb24SBarry Smith     /* create work vectors for each split */
3681b9fc7fcSBarry Smith     if (!jac->x) {
36997bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
3705a9f2f41SSatish Balay       ilink = jac->head;
37197bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
372906ed7ccSBarry Smith         Vec *vl,*vr;
373906ed7ccSBarry Smith 
374906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
375906ed7ccSBarry Smith         ilink->x  = *vr;
376906ed7ccSBarry Smith         ilink->y  = *vl;
377906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
378906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
3795a9f2f41SSatish Balay         jac->x[i] = ilink->x;
3805a9f2f41SSatish Balay         jac->y[i] = ilink->y;
3815a9f2f41SSatish Balay         ilink     = ilink->next;
38297bbdb24SBarry Smith       }
3833b224e63SBarry Smith     }
3843b224e63SBarry Smith   }
3853b224e63SBarry Smith 
3863b224e63SBarry Smith 
3873b224e63SBarry Smith   if (1) {
3883b224e63SBarry Smith     Vec xtmp;
3893b224e63SBarry Smith 
39079416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
3911b9fc7fcSBarry Smith 
3925a9f2f41SSatish Balay     ilink = jac->head;
3931b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
3941b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
395704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
3965a9f2f41SSatish Balay       ilink = ilink->next;
3971b9fc7fcSBarry Smith     }
3981b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
3991b9fc7fcSBarry Smith   }
4000971522cSBarry Smith   PetscFunctionReturn(0);
4010971522cSBarry Smith }
4020971522cSBarry Smith 
4035a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
404ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
405ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
4065a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
407ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
408ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
40979416396SBarry Smith 
4100971522cSBarry Smith #undef __FUNCT__
4113b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
4123b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
4133b224e63SBarry Smith {
4143b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4153b224e63SBarry Smith   PetscErrorCode    ierr;
4163b224e63SBarry Smith   KSP               ksp;
4173b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
4183b224e63SBarry Smith 
4193b224e63SBarry Smith   PetscFunctionBegin;
4203b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
4213b224e63SBarry Smith 
4223b224e63SBarry Smith   ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4233b224e63SBarry Smith   ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4243b224e63SBarry Smith   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
4253b224e63SBarry Smith   ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
4263b224e63SBarry Smith   ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
4273b224e63SBarry Smith   ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4283b224e63SBarry Smith   ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4293b224e63SBarry Smith 
4303b224e63SBarry Smith   ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
4313b224e63SBarry Smith   ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4323b224e63SBarry Smith   ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4333b224e63SBarry Smith 
4343b224e63SBarry Smith   ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
4353b224e63SBarry Smith   ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
4363b224e63SBarry Smith   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
4373b224e63SBarry Smith   ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4383b224e63SBarry Smith   ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4393b224e63SBarry Smith 
4403b224e63SBarry Smith   PetscFunctionReturn(0);
4413b224e63SBarry Smith }
4423b224e63SBarry Smith 
4433b224e63SBarry Smith #undef __FUNCT__
4440971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
4450971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
4460971522cSBarry Smith {
4470971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4480971522cSBarry Smith   PetscErrorCode    ierr;
4495a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
4503e197d65SBarry Smith   PetscInt          bs,cnt;
4510971522cSBarry Smith 
4520971522cSBarry Smith   PetscFunctionBegin;
45351f519a2SBarry Smith   CHKMEMQ;
454e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
45551f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
45651f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
45751f519a2SBarry Smith 
45879416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
4591b9fc7fcSBarry Smith     if (jac->defaultsplit) {
4600971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
4615a9f2f41SSatish Balay       while (ilink) {
4625a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
4635a9f2f41SSatish Balay         ilink = ilink->next;
4640971522cSBarry Smith       }
4650971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
4661b9fc7fcSBarry Smith     } else {
467efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
4685a9f2f41SSatish Balay       while (ilink) {
4695a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
4705a9f2f41SSatish Balay         ilink = ilink->next;
4711b9fc7fcSBarry Smith       }
4721b9fc7fcSBarry Smith     }
47316913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
47479416396SBarry Smith     if (!jac->w1) {
47579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
47679416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
47779416396SBarry Smith     }
478efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
4795a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
4803e197d65SBarry Smith     cnt = 1;
4815a9f2f41SSatish Balay     while (ilink->next) {
4825a9f2f41SSatish Balay       ilink = ilink->next;
4833e197d65SBarry Smith       if (jac->Afield) {
4843e197d65SBarry Smith         /* compute the residual only over the part of the vector needed */
4853e197d65SBarry Smith         ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
4863e197d65SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
4873e197d65SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4883e197d65SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4893e197d65SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
4903e197d65SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4913e197d65SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4923e197d65SBarry Smith       } else {
4933e197d65SBarry Smith         /* compute the residual over the entire vector */
4941ab39975SBarry Smith 	ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
495efb30889SBarry Smith 	ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
4965a9f2f41SSatish Balay 	ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
49779416396SBarry Smith       }
4983e197d65SBarry Smith     }
49951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
50011755939SBarry Smith       cnt -= 2;
50151f519a2SBarry Smith       while (ilink->previous) {
50251f519a2SBarry Smith         ilink = ilink->previous;
50311755939SBarry Smith         if (jac->Afield) {
50411755939SBarry Smith 	  /* compute the residual only over the part of the vector needed */
50511755939SBarry Smith 	  ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
50611755939SBarry Smith 	  ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
50711755939SBarry Smith 	  ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
50811755939SBarry Smith 	  ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
50911755939SBarry Smith 	  ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
51011755939SBarry Smith 	  ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
51111755939SBarry Smith 	  ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
51211755939SBarry Smith         } else {
5131ab39975SBarry Smith 	  ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
51451f519a2SBarry Smith 	  ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
51551f519a2SBarry Smith 	  ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
51679416396SBarry Smith         }
51751f519a2SBarry Smith       }
51811755939SBarry Smith     }
51916913363SBarry Smith   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
52051f519a2SBarry Smith   CHKMEMQ;
5210971522cSBarry Smith   PetscFunctionReturn(0);
5220971522cSBarry Smith }
5230971522cSBarry Smith 
524421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
525ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
526ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
527421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
528ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
529ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
530421e10b8SBarry Smith 
531421e10b8SBarry Smith #undef __FUNCT__
532421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
533421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
534421e10b8SBarry Smith {
535421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
536421e10b8SBarry Smith   PetscErrorCode    ierr;
537421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
538421e10b8SBarry Smith   PetscInt          bs;
539421e10b8SBarry Smith 
540421e10b8SBarry Smith   PetscFunctionBegin;
541421e10b8SBarry Smith   CHKMEMQ;
542421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
543421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
544421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
545421e10b8SBarry Smith 
546421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
547421e10b8SBarry Smith     if (jac->defaultsplit) {
548421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
549421e10b8SBarry Smith       while (ilink) {
550421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
551421e10b8SBarry Smith 	ilink = ilink->next;
552421e10b8SBarry Smith       }
553421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
554421e10b8SBarry Smith     } else {
555421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
556421e10b8SBarry Smith       while (ilink) {
557421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
558421e10b8SBarry Smith 	ilink = ilink->next;
559421e10b8SBarry Smith       }
560421e10b8SBarry Smith     }
561421e10b8SBarry Smith   } else {
562421e10b8SBarry Smith     if (!jac->w1) {
563421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
564421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
565421e10b8SBarry Smith     }
566421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
567421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
568421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
569421e10b8SBarry Smith       while (ilink->next) {
570421e10b8SBarry Smith         ilink = ilink->next;
5719989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
572421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
573421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
574421e10b8SBarry Smith       }
575421e10b8SBarry Smith       while (ilink->previous) {
576421e10b8SBarry Smith         ilink = ilink->previous;
5779989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
578421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
579421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
580421e10b8SBarry Smith       }
581421e10b8SBarry Smith     } else {
582421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
583421e10b8SBarry Smith 	ilink = ilink->next;
584421e10b8SBarry Smith       }
585421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
586421e10b8SBarry Smith       while (ilink->previous) {
587421e10b8SBarry Smith 	ilink = ilink->previous;
5889989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
589421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
590421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
591421e10b8SBarry Smith       }
592421e10b8SBarry Smith     }
593421e10b8SBarry Smith   }
594421e10b8SBarry Smith   CHKMEMQ;
595421e10b8SBarry Smith   PetscFunctionReturn(0);
596421e10b8SBarry Smith }
597421e10b8SBarry Smith 
5980971522cSBarry Smith #undef __FUNCT__
5990971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
6000971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
6010971522cSBarry Smith {
6020971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6030971522cSBarry Smith   PetscErrorCode    ierr;
6045a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
6050971522cSBarry Smith 
6060971522cSBarry Smith   PetscFunctionBegin;
6075a9f2f41SSatish Balay   while (ilink) {
6085a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
6095a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
6105a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
6115a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
612704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
613704ba839SBarry Smith     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
6145a9f2f41SSatish Balay     next = ilink->next;
615704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
616704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
6175a9f2f41SSatish Balay     ilink = next;
6180971522cSBarry Smith   }
61905b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
62097bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
62168dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
62279416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
62379416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
6243b224e63SBarry Smith   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
6253b224e63SBarry Smith   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
6263b224e63SBarry Smith   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
6270971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
6280971522cSBarry Smith   PetscFunctionReturn(0);
6290971522cSBarry Smith }
6300971522cSBarry Smith 
6310971522cSBarry Smith #undef __FUNCT__
6320971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
6330971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
6340971522cSBarry Smith {
6351b9fc7fcSBarry Smith   PetscErrorCode  ierr;
63651f519a2SBarry Smith   PetscInt        i = 0,nfields,*fields,bs;
637e69d4d44SBarry Smith   PetscTruth      flg,set;
6381b9fc7fcSBarry Smith   char            optionname[128];
6399dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
6403b224e63SBarry Smith   PCCompositeType ctype;
6411b9fc7fcSBarry Smith 
6420971522cSBarry Smith   PetscFunctionBegin;
6431b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
64451f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
64551f519a2SBarry Smith   if (flg) {
64651f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
64751f519a2SBarry Smith   }
648704ba839SBarry Smith 
6493b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
6503b224e63SBarry Smith   if (flg) {
6513b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
6523b224e63SBarry Smith   }
653d32f9abdSBarry Smith 
654d32f9abdSBarry Smith   if (jac->bs > 0) {
655d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
656d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
65751f519a2SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
6581b9fc7fcSBarry Smith     while (PETSC_TRUE) {
65913f74950SBarry Smith       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
66051f519a2SBarry Smith       nfields = jac->bs;
6611b9fc7fcSBarry Smith       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
6621b9fc7fcSBarry Smith       if (!flg) break;
6631b9fc7fcSBarry Smith       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
6641b9fc7fcSBarry Smith       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
6651b9fc7fcSBarry Smith     }
66651f519a2SBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
667d32f9abdSBarry Smith   }
668e69d4d44SBarry Smith   ierr = PetscOptionsTruth("-pc_fieldsplit_schur_precondition","Build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",jac->schurpre,&flg,&set);CHKERRQ(ierr);
669e69d4d44SBarry Smith   if (set) {
670e69d4d44SBarry Smith     ierr = PCFieldSplitSchurPrecondition(pc,flg);CHKERRQ(ierr);
671e69d4d44SBarry Smith   }
6721b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6730971522cSBarry Smith   PetscFunctionReturn(0);
6740971522cSBarry Smith }
6750971522cSBarry Smith 
6760971522cSBarry Smith /*------------------------------------------------------------------------------------*/
6770971522cSBarry Smith 
6780971522cSBarry Smith EXTERN_C_BEGIN
6790971522cSBarry Smith #undef __FUNCT__
6800971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
681dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
6820971522cSBarry Smith {
68397bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6840971522cSBarry Smith   PetscErrorCode    ierr;
6855a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
68669a612a9SBarry Smith   char              prefix[128];
68751f519a2SBarry Smith   PetscInt          i;
6880971522cSBarry Smith 
6890971522cSBarry Smith   PetscFunctionBegin;
6900971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
69151f519a2SBarry Smith   for (i=0; i<n; i++) {
69251f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
69351f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
69451f519a2SBarry Smith   }
695704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
696704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
6975a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
6985a9f2f41SSatish Balay   ilink->nfields = n;
6995a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
7007adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
7011cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7025a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
70369a612a9SBarry Smith 
7047adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
7057adad957SLisandro Dalcin     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
70669a612a9SBarry Smith   } else {
70713f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
70869a612a9SBarry Smith   }
7095a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
7100971522cSBarry Smith 
7110971522cSBarry Smith   if (!next) {
7125a9f2f41SSatish Balay     jac->head       = ilink;
71351f519a2SBarry Smith     ilink->previous = PETSC_NULL;
7140971522cSBarry Smith   } else {
7150971522cSBarry Smith     while (next->next) {
7160971522cSBarry Smith       next = next->next;
7170971522cSBarry Smith     }
7185a9f2f41SSatish Balay     next->next      = ilink;
71951f519a2SBarry Smith     ilink->previous = next;
7200971522cSBarry Smith   }
7210971522cSBarry Smith   jac->nsplits++;
7220971522cSBarry Smith   PetscFunctionReturn(0);
7230971522cSBarry Smith }
7240971522cSBarry Smith EXTERN_C_END
7250971522cSBarry Smith 
726e69d4d44SBarry Smith EXTERN_C_BEGIN
727e69d4d44SBarry Smith #undef __FUNCT__
728e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
729e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
730e69d4d44SBarry Smith {
731e69d4d44SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
732e69d4d44SBarry Smith   PetscErrorCode    ierr;
733e69d4d44SBarry Smith   PetscInt          cnt = 0;
734e69d4d44SBarry Smith   PC_FieldSplitLink ilink = jac->head;
735e69d4d44SBarry Smith 
736e69d4d44SBarry Smith   PetscFunctionBegin;
737e69d4d44SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
738e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
739e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
740e69d4d44SBarry Smith   PetscFunctionReturn(0);
741e69d4d44SBarry Smith }
742e69d4d44SBarry Smith EXTERN_C_END
7430971522cSBarry Smith 
7440971522cSBarry Smith EXTERN_C_BEGIN
7450971522cSBarry Smith #undef __FUNCT__
74669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
747dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
7480971522cSBarry Smith {
7490971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7500971522cSBarry Smith   PetscErrorCode    ierr;
7510971522cSBarry Smith   PetscInt          cnt = 0;
7525a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
7530971522cSBarry Smith 
7540971522cSBarry Smith   PetscFunctionBegin;
75569a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
7565a9f2f41SSatish Balay   while (ilink) {
7575a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
7585a9f2f41SSatish Balay     ilink = ilink->next;
7590971522cSBarry Smith   }
7600971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
7610971522cSBarry Smith   *n = jac->nsplits;
7620971522cSBarry Smith   PetscFunctionReturn(0);
7630971522cSBarry Smith }
7640971522cSBarry Smith EXTERN_C_END
7650971522cSBarry Smith 
766704ba839SBarry Smith EXTERN_C_BEGIN
767704ba839SBarry Smith #undef __FUNCT__
768704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
769704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
770704ba839SBarry Smith {
771704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
772704ba839SBarry Smith   PetscErrorCode    ierr;
773704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
774704ba839SBarry Smith   char              prefix[128];
775704ba839SBarry Smith 
776704ba839SBarry Smith   PetscFunctionBegin;
77716913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
7781ab39975SBarry Smith   ilink->is      = is;
7791ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
780704ba839SBarry Smith   ilink->next    = PETSC_NULL;
781704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
7821cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
783704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
784704ba839SBarry Smith 
785704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
786704ba839SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
787704ba839SBarry Smith   } else {
788704ba839SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
789704ba839SBarry Smith   }
790704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
791704ba839SBarry Smith 
792704ba839SBarry Smith   if (!next) {
793704ba839SBarry Smith     jac->head       = ilink;
794704ba839SBarry Smith     ilink->previous = PETSC_NULL;
795704ba839SBarry Smith   } else {
796704ba839SBarry Smith     while (next->next) {
797704ba839SBarry Smith       next = next->next;
798704ba839SBarry Smith     }
799704ba839SBarry Smith     next->next      = ilink;
800704ba839SBarry Smith     ilink->previous = next;
801704ba839SBarry Smith   }
802704ba839SBarry Smith   jac->nsplits++;
803704ba839SBarry Smith 
804704ba839SBarry Smith   PetscFunctionReturn(0);
805704ba839SBarry Smith }
806704ba839SBarry Smith EXTERN_C_END
807704ba839SBarry Smith 
8080971522cSBarry Smith #undef __FUNCT__
8090971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
8100971522cSBarry Smith /*@
8110971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
8120971522cSBarry Smith 
8130971522cSBarry Smith     Collective on PC
8140971522cSBarry Smith 
8150971522cSBarry Smith     Input Parameters:
8160971522cSBarry Smith +   pc  - the preconditioner context
8170971522cSBarry Smith .   n - the number of fields in this split
8180971522cSBarry Smith .   fields - the fields in this split
8190971522cSBarry Smith 
8200971522cSBarry Smith     Level: intermediate
8210971522cSBarry Smith 
822d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
823d32f9abdSBarry Smith 
824d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
825d32f9abdSBarry 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
826d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
827d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
828d32f9abdSBarry Smith 
829d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
8300971522cSBarry Smith 
8310971522cSBarry Smith @*/
832dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
8330971522cSBarry Smith {
8340971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
8350971522cSBarry Smith 
8360971522cSBarry Smith   PetscFunctionBegin;
8370971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8380971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
8390971522cSBarry Smith   if (f) {
8400971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
8410971522cSBarry Smith   }
8420971522cSBarry Smith   PetscFunctionReturn(0);
8430971522cSBarry Smith }
8440971522cSBarry Smith 
8450971522cSBarry Smith #undef __FUNCT__
846704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
847704ba839SBarry Smith /*@
848704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
849704ba839SBarry Smith 
850704ba839SBarry Smith     Collective on PC
851704ba839SBarry Smith 
852704ba839SBarry Smith     Input Parameters:
853704ba839SBarry Smith +   pc  - the preconditioner context
854704ba839SBarry Smith .   is - the index set that defines the vector elements in this field
855704ba839SBarry Smith 
856d32f9abdSBarry Smith 
857d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
858d32f9abdSBarry Smith 
859704ba839SBarry Smith     Level: intermediate
860704ba839SBarry Smith 
861704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
862704ba839SBarry Smith 
863704ba839SBarry Smith @*/
864704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
865704ba839SBarry Smith {
866704ba839SBarry Smith   PetscErrorCode ierr,(*f)(PC,IS);
867704ba839SBarry Smith 
868704ba839SBarry Smith   PetscFunctionBegin;
869704ba839SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
870704ba839SBarry Smith   PetscValidHeaderSpecific(is,IS_COOKIE,1);
871704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
872704ba839SBarry Smith   if (f) {
873704ba839SBarry Smith     ierr = (*f)(pc,is);CHKERRQ(ierr);
874704ba839SBarry Smith   }
875704ba839SBarry Smith   PetscFunctionReturn(0);
876704ba839SBarry Smith }
877704ba839SBarry Smith 
878704ba839SBarry Smith #undef __FUNCT__
87951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
88051f519a2SBarry Smith /*@
88151f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
88251f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
88351f519a2SBarry Smith 
88451f519a2SBarry Smith     Collective on PC
88551f519a2SBarry Smith 
88651f519a2SBarry Smith     Input Parameters:
88751f519a2SBarry Smith +   pc  - the preconditioner context
88851f519a2SBarry Smith -   bs - the block size
88951f519a2SBarry Smith 
89051f519a2SBarry Smith     Level: intermediate
89151f519a2SBarry Smith 
89251f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
89351f519a2SBarry Smith 
89451f519a2SBarry Smith @*/
89551f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
89651f519a2SBarry Smith {
89751f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
89851f519a2SBarry Smith 
89951f519a2SBarry Smith   PetscFunctionBegin;
90051f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
90151f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
90251f519a2SBarry Smith   if (f) {
90351f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
90451f519a2SBarry Smith   }
90551f519a2SBarry Smith   PetscFunctionReturn(0);
90651f519a2SBarry Smith }
90751f519a2SBarry Smith 
90851f519a2SBarry Smith #undef __FUNCT__
90969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
9100971522cSBarry Smith /*@C
91169a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
9120971522cSBarry Smith 
91369a612a9SBarry Smith    Collective on KSP
9140971522cSBarry Smith 
9150971522cSBarry Smith    Input Parameter:
9160971522cSBarry Smith .  pc - the preconditioner context
9170971522cSBarry Smith 
9180971522cSBarry Smith    Output Parameters:
9190971522cSBarry Smith +  n - the number of split
92069a612a9SBarry Smith -  pc - the array of KSP contexts
9210971522cSBarry Smith 
9220971522cSBarry Smith    Note:
923d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
924d32f9abdSBarry Smith    (not the KSP just the array that contains them).
9250971522cSBarry Smith 
92669a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
9270971522cSBarry Smith 
9280971522cSBarry Smith    Level: advanced
9290971522cSBarry Smith 
9300971522cSBarry Smith .seealso: PCFIELDSPLIT
9310971522cSBarry Smith @*/
932dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
9330971522cSBarry Smith {
93469a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
9350971522cSBarry Smith 
9360971522cSBarry Smith   PetscFunctionBegin;
9370971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
9380971522cSBarry Smith   PetscValidIntPointer(n,2);
93969a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
9400971522cSBarry Smith   if (f) {
94169a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
9420971522cSBarry Smith   } else {
94369a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
9440971522cSBarry Smith   }
9450971522cSBarry Smith   PetscFunctionReturn(0);
9460971522cSBarry Smith }
9470971522cSBarry Smith 
948e69d4d44SBarry Smith #undef __FUNCT__
949e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
950e69d4d44SBarry Smith /*@
951e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
952e69d4d44SBarry Smith       D matrix. Otherwise no preconditioner is used.
953e69d4d44SBarry Smith 
954e69d4d44SBarry Smith     Collective on PC
955e69d4d44SBarry Smith 
956e69d4d44SBarry Smith     Input Parameters:
957e69d4d44SBarry Smith +   pc  - the preconditioner context
958e69d4d44SBarry Smith -   flg - build the preconditioner
959e69d4d44SBarry Smith 
960e69d4d44SBarry Smith     Options Database:
961e69d4d44SBarry Smith .     -pc_fieldsplit_schur_precondition <true,false> default is true
962e69d4d44SBarry Smith 
963e69d4d44SBarry Smith     Level: intermediate
964e69d4d44SBarry Smith 
965e69d4d44SBarry Smith     Notes: What should we do if the user wants to provide a different matrix (like a mass matrix) to use to build the preconditioner
966e69d4d44SBarry Smith 
967e69d4d44SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFIELDSPLIT
968e69d4d44SBarry Smith 
969e69d4d44SBarry Smith @*/
970e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PetscTruth flg)
971e69d4d44SBarry Smith {
972e69d4d44SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
973e69d4d44SBarry Smith 
974e69d4d44SBarry Smith   PetscFunctionBegin;
975e69d4d44SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
976e69d4d44SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr);
977e69d4d44SBarry Smith   if (f) {
978e69d4d44SBarry Smith     ierr = (*f)(pc,flg);CHKERRQ(ierr);
979e69d4d44SBarry Smith   }
980e69d4d44SBarry Smith   PetscFunctionReturn(0);
981e69d4d44SBarry Smith }
982e69d4d44SBarry Smith 
983e69d4d44SBarry Smith EXTERN_C_BEGIN
984e69d4d44SBarry Smith #undef __FUNCT__
985e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
986e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PetscTruth flg)
987e69d4d44SBarry Smith {
988e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
989e69d4d44SBarry Smith 
990e69d4d44SBarry Smith   PetscFunctionBegin;
991e69d4d44SBarry Smith   jac->schurpre = flg;
992e69d4d44SBarry Smith   PetscFunctionReturn(0);
993e69d4d44SBarry Smith }
994e69d4d44SBarry Smith EXTERN_C_END
995e69d4d44SBarry Smith 
99679416396SBarry Smith EXTERN_C_BEGIN
99779416396SBarry Smith #undef __FUNCT__
99879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
999dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
100079416396SBarry Smith {
100179416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1002e69d4d44SBarry Smith   PetscErrorCode ierr;
100379416396SBarry Smith 
100479416396SBarry Smith   PetscFunctionBegin;
100579416396SBarry Smith   jac->type = type;
10063b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
10073b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
10083b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1009e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1010e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1011e69d4d44SBarry Smith 
10123b224e63SBarry Smith   } else {
10133b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
10143b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1015e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1016*9e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
10173b224e63SBarry Smith   }
101879416396SBarry Smith   PetscFunctionReturn(0);
101979416396SBarry Smith }
102079416396SBarry Smith EXTERN_C_END
102179416396SBarry Smith 
102251f519a2SBarry Smith EXTERN_C_BEGIN
102351f519a2SBarry Smith #undef __FUNCT__
102451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
102551f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
102651f519a2SBarry Smith {
102751f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
102851f519a2SBarry Smith 
102951f519a2SBarry Smith   PetscFunctionBegin;
103051f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
103151f519a2SBarry 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);
103251f519a2SBarry Smith   jac->bs = bs;
103351f519a2SBarry Smith   PetscFunctionReturn(0);
103451f519a2SBarry Smith }
103551f519a2SBarry Smith EXTERN_C_END
103651f519a2SBarry Smith 
103779416396SBarry Smith #undef __FUNCT__
103879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1039bc08b0f1SBarry Smith /*@
104079416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
104179416396SBarry Smith 
104279416396SBarry Smith    Collective on PC
104379416396SBarry Smith 
104479416396SBarry Smith    Input Parameter:
104579416396SBarry Smith .  pc - the preconditioner context
1046ce9499c7SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
104779416396SBarry Smith 
104879416396SBarry Smith    Options Database Key:
1049ce9499c7SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
105079416396SBarry Smith 
105179416396SBarry Smith    Level: Developer
105279416396SBarry Smith 
105379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
105479416396SBarry Smith 
105579416396SBarry Smith .seealso: PCCompositeSetType()
105679416396SBarry Smith 
105779416396SBarry Smith @*/
1058dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
105979416396SBarry Smith {
106079416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
106179416396SBarry Smith 
106279416396SBarry Smith   PetscFunctionBegin;
106379416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
106479416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
106579416396SBarry Smith   if (f) {
106679416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
106779416396SBarry Smith   }
106879416396SBarry Smith   PetscFunctionReturn(0);
106979416396SBarry Smith }
107079416396SBarry Smith 
10710971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
10720971522cSBarry Smith /*MC
1073a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
10740971522cSBarry Smith                   fields or groups of fields
10750971522cSBarry Smith 
10760971522cSBarry Smith 
10770971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
1078d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
10790971522cSBarry Smith 
1080a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
108169a612a9SBarry Smith          and set the options directly on the resulting KSP object
10820971522cSBarry Smith 
10830971522cSBarry Smith    Level: intermediate
10840971522cSBarry Smith 
108579416396SBarry Smith    Options Database Keys:
10862e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
10872e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
10882e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
108951f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
1090e69d4d44SBarry Smith .    -pc_splitfield_type <additive,multiplicative,schur,symmetric_multiplicative>
1091e69d4d44SBarry Smith .    -pc_fieldsplit_schur_precondition <true,false> default is true
109279416396SBarry Smith 
10933b224e63SBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -schurAblock_ and -schur,
10943b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
10953b224e63SBarry Smith 
10963b224e63SBarry Smith 
1097d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1098d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1099d32f9abdSBarry Smith 
1100d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1101d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1102d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1103d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1104d32f9abdSBarry Smith 
1105d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
1106d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
1107d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
1108d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1109d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
1110d32f9abdSBarry Smith 
1111e69d4d44SBarry Smith      For the Schur complement preconditioner if J = ( A B )
1112e69d4d44SBarry Smith                                                     ( C D )
1113e69d4d44SBarry Smith      the preconditioner is
1114e69d4d44SBarry Smith               (I   -B inv(A)) ( inv(A)   0    ) (I         0  )
1115e69d4d44SBarry Smith               (0    I       ) (   0    inv(S) ) (-C inv(A) I  )
1116e69d4d44SBarry Smith      where the action of inv(A) is applied using the KSP solver with prefix -schurAblock_. The action of
1117e69d4d44SBarry Smith      inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is
1118e69d4d44SBarry Smith      0 it returns the KSP associated with -schurAblock_ while field number 1 gives -schur_ KSP. By default
1119e69d4d44SBarry Smith      D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1120e69d4d44SBarry Smith      option.
1121e69d4d44SBarry Smith 
11220971522cSBarry Smith    Concepts: physics based preconditioners
11230971522cSBarry Smith 
11240971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1125e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
11260971522cSBarry Smith M*/
11270971522cSBarry Smith 
11280971522cSBarry Smith EXTERN_C_BEGIN
11290971522cSBarry Smith #undef __FUNCT__
11300971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
1131dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
11320971522cSBarry Smith {
11330971522cSBarry Smith   PetscErrorCode ierr;
11340971522cSBarry Smith   PC_FieldSplit  *jac;
11350971522cSBarry Smith 
11360971522cSBarry Smith   PetscFunctionBegin;
113738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
11380971522cSBarry Smith   jac->bs        = -1;
11390971522cSBarry Smith   jac->nsplits   = 0;
11403e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1141e69d4d44SBarry Smith   jac->schurpre  = PETSC_TRUE;
114251f519a2SBarry Smith 
11430971522cSBarry Smith   pc->data     = (void*)jac;
11440971522cSBarry Smith 
11450971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1146421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
11470971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
11480971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
11490971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
11500971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
11510971522cSBarry Smith   pc->ops->applyrichardson   = 0;
11520971522cSBarry Smith 
115369a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
115469a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
11550971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
11560971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1157704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1158704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
115979416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
116079416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
116151f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
116251f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
11630971522cSBarry Smith   PetscFunctionReturn(0);
11640971522cSBarry Smith }
11650971522cSBarry Smith EXTERN_C_END
11660971522cSBarry Smith 
11670971522cSBarry Smith 
1168