xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision db4c96c1a48fe024f59866c99a783f61073b2340)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
36356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
40971522cSBarry Smith 
5084e4875SJed Brown const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6084e4875SJed Brown 
70971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
80971522cSBarry Smith struct _PC_FieldSplitLink {
969a612a9SBarry Smith   KSP               ksp;
100971522cSBarry Smith   Vec               x,y;
11*db4c96c1SJed Brown   char              *splitname;
120971522cSBarry Smith   PetscInt          nfields;
130971522cSBarry Smith   PetscInt          *fields;
141b9fc7fcSBarry Smith   VecScatter        sctx;
154aa3045dSJed Brown   IS                is;
1651f519a2SBarry Smith   PC_FieldSplitLink next,previous;
170971522cSBarry Smith };
180971522cSBarry Smith 
190971522cSBarry Smith typedef struct {
2068dd23aaSBarry Smith   PCCompositeType   type;
21a4efd8eaSMatthew Knepley   PetscTruth        defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
22519d70e2SJed Brown   PetscTruth        realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
2330ad9308SMatthew Knepley   PetscInt          bs;           /* Block size for IS and Mat structures */
2430ad9308SMatthew Knepley   PetscInt          nsplits;      /* Number of field divisions defined */
2579416396SBarry Smith   Vec               *x,*y,w1,w2;
26519d70e2SJed Brown   Mat               *mat;         /* The diagonal block for each split */
27519d70e2SJed Brown   Mat               *pmat;        /* The preconditioning diagonal block for each split */
2830ad9308SMatthew Knepley   Mat               *Afield;      /* The rows of the matrix associated with each split */
29704ba839SBarry Smith   PetscTruth        issetup;
3030ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3130ad9308SMatthew Knepley   Mat               B;            /* The (0,1) block */
3230ad9308SMatthew Knepley   Mat               C;            /* The (1,0) block */
3330ad9308SMatthew Knepley   Mat               schur;        /* The Schur complement S = D - C A^{-1} B */
34084e4875SJed Brown   Mat               schur_user;   /* User-provided preconditioning matrix for the Schur complement */
35084e4875SJed Brown   PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
3630ad9308SMatthew Knepley   KSP               kspschur;     /* The solver for S */
3797bbdb24SBarry Smith   PC_FieldSplitLink head;
380971522cSBarry Smith } PC_FieldSplit;
390971522cSBarry Smith 
4016913363SBarry Smith /*
4116913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
4216913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
4316913363SBarry Smith    PC you could change this.
4416913363SBarry Smith */
45084e4875SJed Brown 
46e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
47084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
48084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
49084e4875SJed Brown {
50084e4875SJed Brown   switch (jac->schurpre) {
51084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
52084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
53084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
54084e4875SJed Brown     default:
55084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
56084e4875SJed Brown   }
57084e4875SJed Brown }
58084e4875SJed Brown 
59084e4875SJed Brown 
600971522cSBarry Smith #undef __FUNCT__
610971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
620971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
630971522cSBarry Smith {
640971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
650971522cSBarry Smith   PetscErrorCode    ierr;
660971522cSBarry Smith   PetscTruth        iascii;
670971522cSBarry Smith   PetscInt          i,j;
685a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
690971522cSBarry Smith 
700971522cSBarry Smith   PetscFunctionBegin;
710971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
720971522cSBarry Smith   if (iascii) {
7351f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
7469a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
750971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
760971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
771ab39975SBarry Smith       if (ilink->fields) {
780971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
7979416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
805a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
8179416396SBarry Smith 	  if (j > 0) {
8279416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
8379416396SBarry Smith 	  }
845a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
850971522cSBarry Smith 	}
860971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
8779416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
881ab39975SBarry Smith       } else {
891ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
901ab39975SBarry Smith       }
915a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
925a9f2f41SSatish Balay       ilink = ilink->next;
930971522cSBarry Smith     }
940971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
950971522cSBarry Smith   } else {
9665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
970971522cSBarry Smith   }
980971522cSBarry Smith   PetscFunctionReturn(0);
990971522cSBarry Smith }
1000971522cSBarry Smith 
1010971522cSBarry Smith #undef __FUNCT__
1023b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1033b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1043b224e63SBarry Smith {
1053b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1063b224e63SBarry Smith   PetscErrorCode    ierr;
1073b224e63SBarry Smith   PetscTruth        iascii;
1083b224e63SBarry Smith   PetscInt          i,j;
1093b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1103b224e63SBarry Smith   KSP               ksp;
1113b224e63SBarry Smith 
1123b224e63SBarry Smith   PetscFunctionBegin;
1133b224e63SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1143b224e63SBarry Smith   if (iascii) {
1153b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr);
1163b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1173b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1183b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1193b224e63SBarry Smith       if (ilink->fields) {
1203b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1213b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1223b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1233b224e63SBarry Smith 	  if (j > 0) {
1243b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1253b224e63SBarry Smith 	  }
1263b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1273b224e63SBarry Smith 	}
1283b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1293b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1303b224e63SBarry Smith       } else {
1313b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1323b224e63SBarry Smith       }
1333b224e63SBarry Smith       ilink = ilink->next;
1343b224e63SBarry Smith     }
1353b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
1363b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
13712cae6f2SJed Brown     if (jac->schur) {
13812cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1393b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
14012cae6f2SJed Brown     } else {
14112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
14212cae6f2SJed Brown     }
1433b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1443b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
1453b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14612cae6f2SJed Brown     if (jac->kspschur) {
1473b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
14812cae6f2SJed Brown     } else {
14912cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15012cae6f2SJed Brown     }
1513b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1523b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1533b224e63SBarry Smith   } else {
15465e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1553b224e63SBarry Smith   }
1563b224e63SBarry Smith   PetscFunctionReturn(0);
1573b224e63SBarry Smith }
1583b224e63SBarry Smith 
1593b224e63SBarry Smith #undef __FUNCT__
16069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
16169a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
1620971522cSBarry Smith {
1630971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
1640971522cSBarry Smith   PetscErrorCode    ierr;
1655a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
166d32f9abdSBarry Smith   PetscInt          i = 0,*ifields,nfields;
167*db4c96c1SJed Brown   PetscTruth        flg = PETSC_FALSE,flg2;
168*db4c96c1SJed Brown   char              optionname[128],splitname[8];
1690971522cSBarry Smith 
1700971522cSBarry Smith   PetscFunctionBegin;
171d32f9abdSBarry Smith   if (!ilink) {
172d32f9abdSBarry Smith 
173521d7252SBarry Smith     if (jac->bs <= 0) {
174704ba839SBarry Smith       if (pc->pmat) {
175521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
176704ba839SBarry Smith       } else {
177704ba839SBarry Smith         jac->bs = 1;
178704ba839SBarry Smith       }
179521d7252SBarry Smith     }
180d32f9abdSBarry Smith 
181d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
182d32f9abdSBarry Smith     if (!flg) {
183d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
184d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
185d32f9abdSBarry Smith       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
186d32f9abdSBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
187d32f9abdSBarry Smith       while (PETSC_TRUE) {
188*db4c96c1SJed Brown         sprintf(splitname,"%d_",(int)i);
189d32f9abdSBarry Smith         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
190d32f9abdSBarry Smith         nfields = jac->bs;
191d32f9abdSBarry Smith         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
192d32f9abdSBarry Smith         if (!flg2) break;
193e32f2f54SBarry Smith         if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
194d32f9abdSBarry Smith         flg = PETSC_FALSE;
195*db4c96c1SJed Brown         ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
196d32f9abdSBarry Smith       }
197d32f9abdSBarry Smith       ierr = PetscFree(ifields);CHKERRQ(ierr);
198d32f9abdSBarry Smith     }
199d32f9abdSBarry Smith 
200d32f9abdSBarry Smith     if (flg) {
201d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
202*db4c96c1SJed Brown       for (i=0; i<jac->bs; i++) {
203*db4c96c1SJed Brown         sprintf(splitname,"%d",(int)i);
204*db4c96c1SJed Brown         ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
20579416396SBarry Smith       }
20697bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
207521d7252SBarry Smith     }
208edf189efSBarry Smith   } else if (jac->nsplits == 1) {
209edf189efSBarry Smith     if (ilink->is) {
210edf189efSBarry Smith       IS       is2;
211edf189efSBarry Smith       PetscInt nmin,nmax;
212edf189efSBarry Smith 
213edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
214edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
215*db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
216edf189efSBarry Smith       ierr = ISDestroy(is2);CHKERRQ(ierr);
217e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
218edf189efSBarry Smith   }
219e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
22069a612a9SBarry Smith   PetscFunctionReturn(0);
22169a612a9SBarry Smith }
22269a612a9SBarry Smith 
22369a612a9SBarry Smith 
22469a612a9SBarry Smith #undef __FUNCT__
22569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
22669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
22769a612a9SBarry Smith {
22869a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
22969a612a9SBarry Smith   PetscErrorCode    ierr;
2305a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
231cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
23269a612a9SBarry Smith   MatStructure      flag = pc->flag;
2336c8605c2SJed Brown   PetscTruth        sorted;
23469a612a9SBarry Smith 
23569a612a9SBarry Smith   PetscFunctionBegin;
23669a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
23797bbdb24SBarry Smith   nsplit = jac->nsplits;
2385a9f2f41SSatish Balay   ilink  = jac->head;
23997bbdb24SBarry Smith 
24097bbdb24SBarry Smith   /* get the matrices for each split */
241704ba839SBarry Smith   if (!jac->issetup) {
2421b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
24397bbdb24SBarry Smith 
244704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
245704ba839SBarry Smith 
246704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
24751f519a2SBarry Smith     bs     = jac->bs;
24897bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
249cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
2501b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
2511b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
2521b9fc7fcSBarry Smith       if (jac->defaultsplit) {
253704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
254704ba839SBarry Smith       } else if (!ilink->is) {
255ccb205f8SBarry Smith         if (ilink->nfields > 1) {
2565a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
2575a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
2581b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
2591b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
2601b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
26197bbdb24SBarry Smith             }
26297bbdb24SBarry Smith           }
263704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
264ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
265ccb205f8SBarry Smith         } else {
266704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
267ccb205f8SBarry Smith         }
2683e197d65SBarry Smith       }
269704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
270e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
271704ba839SBarry Smith       ilink = ilink->next;
2721b9fc7fcSBarry Smith     }
2731b9fc7fcSBarry Smith   }
2741b9fc7fcSBarry Smith 
275704ba839SBarry Smith   ilink  = jac->head;
27697bbdb24SBarry Smith   if (!jac->pmat) {
277cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
278cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
2794aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
280704ba839SBarry Smith       ilink = ilink->next;
281cf502942SBarry Smith     }
28297bbdb24SBarry Smith   } else {
283cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
2844aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
285704ba839SBarry Smith       ilink = ilink->next;
286cf502942SBarry Smith     }
28797bbdb24SBarry Smith   }
288519d70e2SJed Brown   if (jac->realdiagonal) {
289519d70e2SJed Brown     ilink = jac->head;
290519d70e2SJed Brown     if (!jac->mat) {
291519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
292519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
293519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
294519d70e2SJed Brown         ilink = ilink->next;
295519d70e2SJed Brown       }
296519d70e2SJed Brown     } else {
297519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
298519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
299519d70e2SJed Brown         ilink = ilink->next;
300519d70e2SJed Brown       }
301519d70e2SJed Brown     }
302519d70e2SJed Brown   } else {
303519d70e2SJed Brown     jac->mat = jac->pmat;
304519d70e2SJed Brown   }
30597bbdb24SBarry Smith 
3066c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
30768dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
30868dd23aaSBarry Smith     ilink  = jac->head;
30968dd23aaSBarry Smith     if (!jac->Afield) {
31068dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
31168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3124aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
31368dd23aaSBarry Smith         ilink = ilink->next;
31468dd23aaSBarry Smith       }
31568dd23aaSBarry Smith     } else {
31668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3174aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
31868dd23aaSBarry Smith         ilink = ilink->next;
31968dd23aaSBarry Smith       }
32068dd23aaSBarry Smith     }
32168dd23aaSBarry Smith   }
32268dd23aaSBarry Smith 
3233b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
3243b224e63SBarry Smith     IS       ccis;
3254aa3045dSJed Brown     PetscInt rstart,rend;
326e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
32768dd23aaSBarry Smith 
328e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
329e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
330e6cab6aaSJed Brown 
3313b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
3323b224e63SBarry Smith     if (jac->schur) {
3333b224e63SBarry Smith       ilink = jac->head;
3344aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3354aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
3363b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3373b224e63SBarry Smith       ilink = ilink->next;
3384aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3394aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
3403b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
341519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
342084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
3433b224e63SBarry Smith 
3443b224e63SBarry Smith      } else {
3451cee3971SBarry Smith       KSP ksp;
3463b224e63SBarry Smith 
3473b224e63SBarry Smith       /* extract the B and C matrices */
3483b224e63SBarry Smith       ilink = jac->head;
3494aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3504aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
3513b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3523b224e63SBarry Smith       ilink = ilink->next;
3534aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3544aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
3553b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
356084e4875SJed Brown       /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */
357519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
3581cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
359e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
3603b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
3613b224e63SBarry Smith 
3623b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
3631cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
364084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
365084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
366084e4875SJed Brown         PC pc;
367cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
368084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
369084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
370e69d4d44SBarry Smith       }
3713b224e63SBarry Smith       ierr  = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr);
372edf189efSBarry Smith       ierr  = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr);
3733b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
3743b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
3753b224e63SBarry Smith 
3763b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
3773b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
3783b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
3793b224e63SBarry Smith       ilink = jac->head;
3803b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
3813b224e63SBarry Smith       ilink = ilink->next;
3823b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
3833b224e63SBarry Smith     }
3843b224e63SBarry Smith   } else {
38597bbdb24SBarry Smith     /* set up the individual PCs */
38697bbdb24SBarry Smith     i    = 0;
3875a9f2f41SSatish Balay     ilink = jac->head;
3885a9f2f41SSatish Balay     while (ilink) {
389519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
3903b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
3915a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
3925a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
39397bbdb24SBarry Smith       i++;
3945a9f2f41SSatish Balay       ilink = ilink->next;
3950971522cSBarry Smith     }
39697bbdb24SBarry Smith 
39797bbdb24SBarry Smith     /* create work vectors for each split */
3981b9fc7fcSBarry Smith     if (!jac->x) {
39997bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
4005a9f2f41SSatish Balay       ilink = jac->head;
40197bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
402906ed7ccSBarry Smith         Vec *vl,*vr;
403906ed7ccSBarry Smith 
404906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
405906ed7ccSBarry Smith         ilink->x  = *vr;
406906ed7ccSBarry Smith         ilink->y  = *vl;
407906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
408906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
4095a9f2f41SSatish Balay         jac->x[i] = ilink->x;
4105a9f2f41SSatish Balay         jac->y[i] = ilink->y;
4115a9f2f41SSatish Balay         ilink     = ilink->next;
41297bbdb24SBarry Smith       }
4133b224e63SBarry Smith     }
4143b224e63SBarry Smith   }
4153b224e63SBarry Smith 
4163b224e63SBarry Smith 
417d0f46423SBarry Smith   if (!jac->head->sctx) {
4183b224e63SBarry Smith     Vec xtmp;
4193b224e63SBarry Smith 
42079416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
4211b9fc7fcSBarry Smith 
4225a9f2f41SSatish Balay     ilink = jac->head;
4231b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
4241b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
425704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
4265a9f2f41SSatish Balay       ilink = ilink->next;
4271b9fc7fcSBarry Smith     }
4281b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
4291b9fc7fcSBarry Smith   }
4300971522cSBarry Smith   PetscFunctionReturn(0);
4310971522cSBarry Smith }
4320971522cSBarry Smith 
4335a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
434ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
435ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
4365a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
437ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
438ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
43979416396SBarry Smith 
4400971522cSBarry Smith #undef __FUNCT__
4413b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
4423b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
4433b224e63SBarry Smith {
4443b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4453b224e63SBarry Smith   PetscErrorCode    ierr;
4463b224e63SBarry Smith   KSP               ksp;
4473b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
4483b224e63SBarry Smith 
4493b224e63SBarry Smith   PetscFunctionBegin;
4503b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
4513b224e63SBarry Smith 
4523b224e63SBarry Smith   ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4533b224e63SBarry Smith   ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4543b224e63SBarry Smith   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
4553b224e63SBarry Smith   ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
4563b224e63SBarry Smith   ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
4573b224e63SBarry Smith   ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4583b224e63SBarry Smith   ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4593b224e63SBarry Smith 
4603b224e63SBarry Smith   ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
4613b224e63SBarry Smith   ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4623b224e63SBarry Smith   ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4633b224e63SBarry Smith 
4643b224e63SBarry Smith   ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
4653b224e63SBarry Smith   ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
4663b224e63SBarry Smith   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
4673b224e63SBarry Smith   ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4683b224e63SBarry Smith   ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4693b224e63SBarry Smith 
4703b224e63SBarry Smith   PetscFunctionReturn(0);
4713b224e63SBarry Smith }
4723b224e63SBarry Smith 
4733b224e63SBarry Smith #undef __FUNCT__
4740971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
4750971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
4760971522cSBarry Smith {
4770971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4780971522cSBarry Smith   PetscErrorCode    ierr;
4795a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
480af93645bSJed Brown   PetscInt          cnt;
4810971522cSBarry Smith 
4820971522cSBarry Smith   PetscFunctionBegin;
48351f519a2SBarry Smith   CHKMEMQ;
48451f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
48551f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
48651f519a2SBarry Smith 
48779416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
4881b9fc7fcSBarry Smith     if (jac->defaultsplit) {
4890971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
4905a9f2f41SSatish Balay       while (ilink) {
4915a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
4925a9f2f41SSatish Balay         ilink = ilink->next;
4930971522cSBarry Smith       }
4940971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
4951b9fc7fcSBarry Smith     } else {
496efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
4975a9f2f41SSatish Balay       while (ilink) {
4985a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
4995a9f2f41SSatish Balay         ilink = ilink->next;
5001b9fc7fcSBarry Smith       }
5011b9fc7fcSBarry Smith     }
50216913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
50379416396SBarry Smith     if (!jac->w1) {
50479416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
50579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
50679416396SBarry Smith     }
507efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
5085a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
5093e197d65SBarry Smith     cnt = 1;
5105a9f2f41SSatish Balay     while (ilink->next) {
5115a9f2f41SSatish Balay       ilink = ilink->next;
5123e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
5133e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
5143e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
5153e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5163e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5173e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
5183e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5193e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5203e197d65SBarry Smith     }
52151f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
52211755939SBarry Smith       cnt -= 2;
52351f519a2SBarry Smith       while (ilink->previous) {
52451f519a2SBarry Smith         ilink = ilink->previous;
52511755939SBarry Smith         /* compute the residual only over the part of the vector needed */
52611755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
52711755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
52811755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
52911755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
53011755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
53111755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
53211755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
53351f519a2SBarry Smith       }
53411755939SBarry Smith     }
53565e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
53651f519a2SBarry Smith   CHKMEMQ;
5370971522cSBarry Smith   PetscFunctionReturn(0);
5380971522cSBarry Smith }
5390971522cSBarry Smith 
540421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
541ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
542ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
543421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
544ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
545ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
546421e10b8SBarry Smith 
547421e10b8SBarry Smith #undef __FUNCT__
548421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
549421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
550421e10b8SBarry Smith {
551421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
552421e10b8SBarry Smith   PetscErrorCode    ierr;
553421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
554421e10b8SBarry Smith 
555421e10b8SBarry Smith   PetscFunctionBegin;
556421e10b8SBarry Smith   CHKMEMQ;
557421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
558421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
559421e10b8SBarry Smith 
560421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
561421e10b8SBarry Smith     if (jac->defaultsplit) {
562421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
563421e10b8SBarry Smith       while (ilink) {
564421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
565421e10b8SBarry Smith 	ilink = ilink->next;
566421e10b8SBarry Smith       }
567421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
568421e10b8SBarry Smith     } else {
569421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
570421e10b8SBarry Smith       while (ilink) {
571421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
572421e10b8SBarry Smith 	ilink = ilink->next;
573421e10b8SBarry Smith       }
574421e10b8SBarry Smith     }
575421e10b8SBarry Smith   } else {
576421e10b8SBarry Smith     if (!jac->w1) {
577421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
578421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
579421e10b8SBarry Smith     }
580421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
581421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
582421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
583421e10b8SBarry Smith       while (ilink->next) {
584421e10b8SBarry Smith         ilink = ilink->next;
5859989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
586421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
587421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
588421e10b8SBarry Smith       }
589421e10b8SBarry Smith       while (ilink->previous) {
590421e10b8SBarry Smith         ilink = ilink->previous;
5919989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
592421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
593421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
594421e10b8SBarry Smith       }
595421e10b8SBarry Smith     } else {
596421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
597421e10b8SBarry Smith 	ilink = ilink->next;
598421e10b8SBarry Smith       }
599421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
600421e10b8SBarry Smith       while (ilink->previous) {
601421e10b8SBarry Smith 	ilink = ilink->previous;
6029989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
603421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
604421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
605421e10b8SBarry Smith       }
606421e10b8SBarry Smith     }
607421e10b8SBarry Smith   }
608421e10b8SBarry Smith   CHKMEMQ;
609421e10b8SBarry Smith   PetscFunctionReturn(0);
610421e10b8SBarry Smith }
611421e10b8SBarry Smith 
6120971522cSBarry Smith #undef __FUNCT__
6130971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
6140971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
6150971522cSBarry Smith {
6160971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6170971522cSBarry Smith   PetscErrorCode    ierr;
6185a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
6190971522cSBarry Smith 
6200971522cSBarry Smith   PetscFunctionBegin;
6215a9f2f41SSatish Balay   while (ilink) {
6225a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
6235a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
6245a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
6255a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
626704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
6275a9f2f41SSatish Balay     next = ilink->next;
628704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
629704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
6305a9f2f41SSatish Balay     ilink = next;
6310971522cSBarry Smith   }
63205b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
633519d70e2SJed Brown   if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);}
63497bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
63568dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
63679416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
63779416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
6383b224e63SBarry Smith   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
639084e4875SJed Brown   if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
640d0f46423SBarry Smith   if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);}
6413b224e63SBarry Smith   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
6423b224e63SBarry Smith   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
6430971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
6440971522cSBarry Smith   PetscFunctionReturn(0);
6450971522cSBarry Smith }
6460971522cSBarry Smith 
6470971522cSBarry Smith #undef __FUNCT__
6480971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
6490971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
6500971522cSBarry Smith {
6511b9fc7fcSBarry Smith   PetscErrorCode  ierr;
65251f519a2SBarry Smith   PetscInt        i = 0,nfields,*fields,bs;
653084e4875SJed Brown   PetscTruth      flg;
654*db4c96c1SJed Brown   char            optionname[128],splitname[8];
6559dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
6563b224e63SBarry Smith   PCCompositeType ctype;
6571b9fc7fcSBarry Smith 
6580971522cSBarry Smith   PetscFunctionBegin;
6591b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
660519d70e2SJed Brown   ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
66151f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
66251f519a2SBarry Smith   if (flg) {
66351f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
66451f519a2SBarry Smith   }
665704ba839SBarry Smith 
6663b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
6673b224e63SBarry Smith   if (flg) {
6683b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
6693b224e63SBarry Smith   }
670d32f9abdSBarry Smith 
671c30613efSMatthew Knepley   /* Only setup fields once */
672c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
673d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
674d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
67551f519a2SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
6761b9fc7fcSBarry Smith     while (PETSC_TRUE) {
677*db4c96c1SJed Brown       sprintf(splitname,"%d",(int)i);
67813f74950SBarry Smith       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
67951f519a2SBarry Smith       nfields = jac->bs;
6801b9fc7fcSBarry Smith       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
6811b9fc7fcSBarry Smith       if (!flg) break;
682e32f2f54SBarry Smith       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
683*db4c96c1SJed Brown       ierr = PCFieldSplitSetFields(pc,splitname,nfields,fields);CHKERRQ(ierr);
6841b9fc7fcSBarry Smith     }
68551f519a2SBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
686d32f9abdSBarry Smith   }
687084e4875SJed Brown   ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,PETSC_NULL);CHKERRQ(ierr);
6881b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6890971522cSBarry Smith   PetscFunctionReturn(0);
6900971522cSBarry Smith }
6910971522cSBarry Smith 
6920971522cSBarry Smith /*------------------------------------------------------------------------------------*/
6930971522cSBarry Smith 
6940971522cSBarry Smith EXTERN_C_BEGIN
6950971522cSBarry Smith #undef __FUNCT__
6960971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
697*db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,PetscInt *fields)
6980971522cSBarry Smith {
69997bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7000971522cSBarry Smith   PetscErrorCode    ierr;
7015a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
70269a612a9SBarry Smith   char              prefix[128];
70351f519a2SBarry Smith   PetscInt          i;
7040971522cSBarry Smith 
7050971522cSBarry Smith   PetscFunctionBegin;
70651f519a2SBarry Smith   for (i=0; i<n; i++) {
707e32f2f54SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
708e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
70951f519a2SBarry Smith   }
710704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
711*db4c96c1SJed Brown   ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
712704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
7135a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
7145a9f2f41SSatish Balay   ilink->nfields = n;
7155a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
7167adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
7171cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7185a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
71969a612a9SBarry Smith 
7207adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
721*db4c96c1SJed Brown     sprintf(prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix,splitname);
72269a612a9SBarry Smith   } else {
723*db4c96c1SJed Brown     sprintf(prefix,"fieldsplit_%s_",splitname);
72469a612a9SBarry Smith   }
7255a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
7260971522cSBarry Smith 
7270971522cSBarry Smith   if (!next) {
7285a9f2f41SSatish Balay     jac->head       = ilink;
72951f519a2SBarry Smith     ilink->previous = PETSC_NULL;
7300971522cSBarry Smith   } else {
7310971522cSBarry Smith     while (next->next) {
7320971522cSBarry Smith       next = next->next;
7330971522cSBarry Smith     }
7345a9f2f41SSatish Balay     next->next      = ilink;
73551f519a2SBarry Smith     ilink->previous = next;
7360971522cSBarry Smith   }
7370971522cSBarry Smith   jac->nsplits++;
7380971522cSBarry Smith   PetscFunctionReturn(0);
7390971522cSBarry Smith }
7400971522cSBarry Smith EXTERN_C_END
7410971522cSBarry Smith 
742e69d4d44SBarry Smith EXTERN_C_BEGIN
743e69d4d44SBarry Smith #undef __FUNCT__
744e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
745e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
746e69d4d44SBarry Smith {
747e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
748e69d4d44SBarry Smith   PetscErrorCode ierr;
749e69d4d44SBarry Smith 
750e69d4d44SBarry Smith   PetscFunctionBegin;
751519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
752e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
753e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
754084e4875SJed Brown   *n = jac->nsplits;
755e69d4d44SBarry Smith   PetscFunctionReturn(0);
756e69d4d44SBarry Smith }
757e69d4d44SBarry Smith EXTERN_C_END
7580971522cSBarry Smith 
7590971522cSBarry Smith EXTERN_C_BEGIN
7600971522cSBarry Smith #undef __FUNCT__
76169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
762dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
7630971522cSBarry Smith {
7640971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7650971522cSBarry Smith   PetscErrorCode    ierr;
7660971522cSBarry Smith   PetscInt          cnt = 0;
7675a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
7680971522cSBarry Smith 
7690971522cSBarry Smith   PetscFunctionBegin;
77069a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
7715a9f2f41SSatish Balay   while (ilink) {
7725a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
7735a9f2f41SSatish Balay     ilink = ilink->next;
7740971522cSBarry Smith   }
775e32f2f54SBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
7760971522cSBarry Smith   *n = jac->nsplits;
7770971522cSBarry Smith   PetscFunctionReturn(0);
7780971522cSBarry Smith }
7790971522cSBarry Smith EXTERN_C_END
7800971522cSBarry Smith 
781704ba839SBarry Smith EXTERN_C_BEGIN
782704ba839SBarry Smith #undef __FUNCT__
783704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
784*db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
785704ba839SBarry Smith {
786704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
787704ba839SBarry Smith   PetscErrorCode    ierr;
788704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
789704ba839SBarry Smith   char              prefix[128];
790704ba839SBarry Smith 
791704ba839SBarry Smith   PetscFunctionBegin;
79216913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
793*db4c96c1SJed Brown   ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
7941ab39975SBarry Smith   ilink->is      = is;
7951ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
796704ba839SBarry Smith   ilink->next    = PETSC_NULL;
797704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
7981cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
799704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
800704ba839SBarry Smith 
801704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
802*db4c96c1SJed Brown     sprintf(prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix,splitname);
803704ba839SBarry Smith   } else {
804*db4c96c1SJed Brown     sprintf(prefix,"fieldsplit_%s_",splitname);
805704ba839SBarry Smith   }
806704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
807704ba839SBarry Smith 
808704ba839SBarry Smith   if (!next) {
809704ba839SBarry Smith     jac->head       = ilink;
810704ba839SBarry Smith     ilink->previous = PETSC_NULL;
811704ba839SBarry Smith   } else {
812704ba839SBarry Smith     while (next->next) {
813704ba839SBarry Smith       next = next->next;
814704ba839SBarry Smith     }
815704ba839SBarry Smith     next->next      = ilink;
816704ba839SBarry Smith     ilink->previous = next;
817704ba839SBarry Smith   }
818704ba839SBarry Smith   jac->nsplits++;
819704ba839SBarry Smith 
820704ba839SBarry Smith   PetscFunctionReturn(0);
821704ba839SBarry Smith }
822704ba839SBarry Smith EXTERN_C_END
823704ba839SBarry Smith 
8240971522cSBarry Smith #undef __FUNCT__
8250971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
8260971522cSBarry Smith /*@
8270971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
8280971522cSBarry Smith 
8290971522cSBarry Smith     Collective on PC
8300971522cSBarry Smith 
8310971522cSBarry Smith     Input Parameters:
8320971522cSBarry Smith +   pc  - the preconditioner context
833*db4c96c1SJed Brown .   splitname - name of this split
8340971522cSBarry Smith .   n - the number of fields in this split
835*db4c96c1SJed Brown -   fields - the fields in this split
8360971522cSBarry Smith 
8370971522cSBarry Smith     Level: intermediate
8380971522cSBarry Smith 
839d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
840d32f9abdSBarry Smith 
841d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
842d32f9abdSBarry 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
843d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
844d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
845d32f9abdSBarry Smith 
846*db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
847*db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
848*db4c96c1SJed Brown 
849d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
8500971522cSBarry Smith 
8510971522cSBarry Smith @*/
852*db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n, PetscInt *fields)
8530971522cSBarry Smith {
854*db4c96c1SJed Brown   PetscErrorCode ierr,(*f)(PC,const char[],PetscInt,PetscInt *);
8550971522cSBarry Smith 
8560971522cSBarry Smith   PetscFunctionBegin;
8570700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
858*db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
859*db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
860*db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
8610971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
8620971522cSBarry Smith   if (f) {
863*db4c96c1SJed Brown     ierr = (*f)(pc,splitname,n,fields);CHKERRQ(ierr);
8640971522cSBarry Smith   }
8650971522cSBarry Smith   PetscFunctionReturn(0);
8660971522cSBarry Smith }
8670971522cSBarry Smith 
8680971522cSBarry Smith #undef __FUNCT__
869704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
870704ba839SBarry Smith /*@
871704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
872704ba839SBarry Smith 
873704ba839SBarry Smith     Collective on PC
874704ba839SBarry Smith 
875704ba839SBarry Smith     Input Parameters:
876704ba839SBarry Smith +   pc  - the preconditioner context
877*db4c96c1SJed Brown .   splitname - name of this split
878*db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
879704ba839SBarry Smith 
880d32f9abdSBarry Smith 
881a6ffb8dbSJed Brown     Notes:
882a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
883a6ffb8dbSJed Brown 
884*db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
885*db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
886d32f9abdSBarry Smith 
887704ba839SBarry Smith     Level: intermediate
888704ba839SBarry Smith 
889704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
890704ba839SBarry Smith 
891704ba839SBarry Smith @*/
892*db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
893704ba839SBarry Smith {
894*db4c96c1SJed Brown   PetscErrorCode ierr,(*f)(PC,const char[],IS);
895704ba839SBarry Smith 
896704ba839SBarry Smith   PetscFunctionBegin;
8970700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
898*db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
899*db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
900704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
901704ba839SBarry Smith   if (f) {
902*db4c96c1SJed Brown     ierr = (*f)(pc,splitname,is);CHKERRQ(ierr);
903704ba839SBarry Smith   }
904704ba839SBarry Smith   PetscFunctionReturn(0);
905704ba839SBarry Smith }
906704ba839SBarry Smith 
907704ba839SBarry Smith #undef __FUNCT__
90851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
90951f519a2SBarry Smith /*@
91051f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
91151f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
91251f519a2SBarry Smith 
91351f519a2SBarry Smith     Collective on PC
91451f519a2SBarry Smith 
91551f519a2SBarry Smith     Input Parameters:
91651f519a2SBarry Smith +   pc  - the preconditioner context
91751f519a2SBarry Smith -   bs - the block size
91851f519a2SBarry Smith 
91951f519a2SBarry Smith     Level: intermediate
92051f519a2SBarry Smith 
92151f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
92251f519a2SBarry Smith 
92351f519a2SBarry Smith @*/
92451f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
92551f519a2SBarry Smith {
92651f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
92751f519a2SBarry Smith 
92851f519a2SBarry Smith   PetscFunctionBegin;
9290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
93051f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
93151f519a2SBarry Smith   if (f) {
93251f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
93351f519a2SBarry Smith   }
93451f519a2SBarry Smith   PetscFunctionReturn(0);
93551f519a2SBarry Smith }
93651f519a2SBarry Smith 
93751f519a2SBarry Smith #undef __FUNCT__
93869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
9390971522cSBarry Smith /*@C
94069a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
9410971522cSBarry Smith 
94269a612a9SBarry Smith    Collective on KSP
9430971522cSBarry Smith 
9440971522cSBarry Smith    Input Parameter:
9450971522cSBarry Smith .  pc - the preconditioner context
9460971522cSBarry Smith 
9470971522cSBarry Smith    Output Parameters:
9480971522cSBarry Smith +  n - the number of split
94969a612a9SBarry Smith -  pc - the array of KSP contexts
9500971522cSBarry Smith 
9510971522cSBarry Smith    Note:
952d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
953d32f9abdSBarry Smith    (not the KSP just the array that contains them).
9540971522cSBarry Smith 
95569a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
9560971522cSBarry Smith 
9570971522cSBarry Smith    Level: advanced
9580971522cSBarry Smith 
9590971522cSBarry Smith .seealso: PCFIELDSPLIT
9600971522cSBarry Smith @*/
961dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
9620971522cSBarry Smith {
96369a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
9640971522cSBarry Smith 
9650971522cSBarry Smith   PetscFunctionBegin;
9660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9670971522cSBarry Smith   PetscValidIntPointer(n,2);
96869a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
9690971522cSBarry Smith   if (f) {
97069a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
971e7e72b3dSBarry Smith   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
9720971522cSBarry Smith   PetscFunctionReturn(0);
9730971522cSBarry Smith }
9740971522cSBarry Smith 
975e69d4d44SBarry Smith #undef __FUNCT__
976e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
977e69d4d44SBarry Smith /*@
978e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
979e69d4d44SBarry Smith       D matrix. Otherwise no preconditioner is used.
980e69d4d44SBarry Smith 
981e69d4d44SBarry Smith     Collective on PC
982e69d4d44SBarry Smith 
983e69d4d44SBarry Smith     Input Parameters:
984e69d4d44SBarry Smith +   pc  - the preconditioner context
985084e4875SJed Brown .   ptype - which matrix to use for preconditioning the Schur complement
986084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
987084e4875SJed Brown 
988084e4875SJed Brown     Notes:
989084e4875SJed Brown     The default is to use the block on the diagonal of the preconditioning matrix.  This is D, in the (1,1) position.
990084e4875SJed Brown     There are currently no preconditioners that work directly with the Schur complement so setting
991084e4875SJed Brown     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
992e69d4d44SBarry Smith 
993e69d4d44SBarry Smith     Options Database:
994084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
995e69d4d44SBarry Smith 
996e69d4d44SBarry Smith     Level: intermediate
997e69d4d44SBarry Smith 
998084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
999e69d4d44SBarry Smith 
1000e69d4d44SBarry Smith @*/
1001084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1002e69d4d44SBarry Smith {
1003084e4875SJed Brown   PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat);
1004e69d4d44SBarry Smith 
1005e69d4d44SBarry Smith   PetscFunctionBegin;
10060700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1007e69d4d44SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr);
1008e69d4d44SBarry Smith   if (f) {
1009084e4875SJed Brown     ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr);
1010e69d4d44SBarry Smith   }
1011e69d4d44SBarry Smith   PetscFunctionReturn(0);
1012e69d4d44SBarry Smith }
1013e69d4d44SBarry Smith 
1014e69d4d44SBarry Smith EXTERN_C_BEGIN
1015e69d4d44SBarry Smith #undef __FUNCT__
1016e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1017084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1018e69d4d44SBarry Smith {
1019e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1020084e4875SJed Brown   PetscErrorCode  ierr;
1021e69d4d44SBarry Smith 
1022e69d4d44SBarry Smith   PetscFunctionBegin;
1023084e4875SJed Brown   jac->schurpre = ptype;
1024084e4875SJed Brown   if (pre) {
1025084e4875SJed Brown     if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
1026084e4875SJed Brown     jac->schur_user = pre;
1027084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1028084e4875SJed Brown   }
1029e69d4d44SBarry Smith   PetscFunctionReturn(0);
1030e69d4d44SBarry Smith }
1031e69d4d44SBarry Smith EXTERN_C_END
1032e69d4d44SBarry Smith 
103330ad9308SMatthew Knepley #undef __FUNCT__
103430ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
103530ad9308SMatthew Knepley /*@C
103630ad9308SMatthew Knepley    PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement
103730ad9308SMatthew Knepley 
103830ad9308SMatthew Knepley    Collective on KSP
103930ad9308SMatthew Knepley 
104030ad9308SMatthew Knepley    Input Parameter:
104130ad9308SMatthew Knepley .  pc - the preconditioner context
104230ad9308SMatthew Knepley 
104330ad9308SMatthew Knepley    Output Parameters:
104430ad9308SMatthew Knepley +  A - the (0,0) block
104530ad9308SMatthew Knepley .  B - the (0,1) block
104630ad9308SMatthew Knepley .  C - the (1,0) block
104730ad9308SMatthew Knepley -  D - the (1,1) block
104830ad9308SMatthew Knepley 
104930ad9308SMatthew Knepley    Level: advanced
105030ad9308SMatthew Knepley 
105130ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
105230ad9308SMatthew Knepley @*/
105330ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D)
105430ad9308SMatthew Knepley {
105530ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
105630ad9308SMatthew Knepley 
105730ad9308SMatthew Knepley   PetscFunctionBegin;
10580700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1059c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
106030ad9308SMatthew Knepley   if (A) *A = jac->pmat[0];
106130ad9308SMatthew Knepley   if (B) *B = jac->B;
106230ad9308SMatthew Knepley   if (C) *C = jac->C;
106330ad9308SMatthew Knepley   if (D) *D = jac->pmat[1];
106430ad9308SMatthew Knepley   PetscFunctionReturn(0);
106530ad9308SMatthew Knepley }
106630ad9308SMatthew Knepley 
106779416396SBarry Smith EXTERN_C_BEGIN
106879416396SBarry Smith #undef __FUNCT__
106979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1070dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
107179416396SBarry Smith {
107279416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1073e69d4d44SBarry Smith   PetscErrorCode ierr;
107479416396SBarry Smith 
107579416396SBarry Smith   PetscFunctionBegin;
107679416396SBarry Smith   jac->type = type;
10773b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
10783b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
10793b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1080e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1081e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1082e69d4d44SBarry Smith 
10833b224e63SBarry Smith   } else {
10843b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
10853b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1086e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
10879e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
10883b224e63SBarry Smith   }
108979416396SBarry Smith   PetscFunctionReturn(0);
109079416396SBarry Smith }
109179416396SBarry Smith EXTERN_C_END
109279416396SBarry Smith 
109351f519a2SBarry Smith EXTERN_C_BEGIN
109451f519a2SBarry Smith #undef __FUNCT__
109551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
109651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
109751f519a2SBarry Smith {
109851f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
109951f519a2SBarry Smith 
110051f519a2SBarry Smith   PetscFunctionBegin;
110165e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
110265e19b50SBarry Smith   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
110351f519a2SBarry Smith   jac->bs = bs;
110451f519a2SBarry Smith   PetscFunctionReturn(0);
110551f519a2SBarry Smith }
110651f519a2SBarry Smith EXTERN_C_END
110751f519a2SBarry Smith 
110879416396SBarry Smith #undef __FUNCT__
110979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1110bc08b0f1SBarry Smith /*@
111179416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
111279416396SBarry Smith 
111379416396SBarry Smith    Collective on PC
111479416396SBarry Smith 
111579416396SBarry Smith    Input Parameter:
111679416396SBarry Smith .  pc - the preconditioner context
111781540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
111879416396SBarry Smith 
111979416396SBarry Smith    Options Database Key:
1120a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
112179416396SBarry Smith 
112279416396SBarry Smith    Level: Developer
112379416396SBarry Smith 
112479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
112579416396SBarry Smith 
112679416396SBarry Smith .seealso: PCCompositeSetType()
112779416396SBarry Smith 
112879416396SBarry Smith @*/
1129dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
113079416396SBarry Smith {
113179416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
113279416396SBarry Smith 
113379416396SBarry Smith   PetscFunctionBegin;
11340700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
113579416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
113679416396SBarry Smith   if (f) {
113779416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
113879416396SBarry Smith   }
113979416396SBarry Smith   PetscFunctionReturn(0);
114079416396SBarry Smith }
114179416396SBarry Smith 
11420971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
11430971522cSBarry Smith /*MC
1144a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
11450971522cSBarry Smith                   fields or groups of fields
11460971522cSBarry Smith 
11470971522cSBarry Smith 
1148edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1149edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
11500971522cSBarry Smith 
1151a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
115269a612a9SBarry Smith          and set the options directly on the resulting KSP object
11530971522cSBarry Smith 
11540971522cSBarry Smith    Level: intermediate
11550971522cSBarry Smith 
115679416396SBarry Smith    Options Database Keys:
115781540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
115881540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
115981540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
116081540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
116181540f2fSBarry Smith .    -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1162e69d4d44SBarry Smith .    -pc_fieldsplit_schur_precondition <true,false> default is true
116379416396SBarry Smith 
1164edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
11653b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
11663b224e63SBarry Smith 
11673b224e63SBarry Smith 
1168d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1169d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1170d32f9abdSBarry Smith 
1171d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1172d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1173d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1174d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1175d32f9abdSBarry Smith 
1176d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
1177d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
1178d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
1179d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1180d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
1181d32f9abdSBarry Smith 
1182e69d4d44SBarry Smith      For the Schur complement preconditioner if J = ( A B )
1183e69d4d44SBarry Smith                                                     ( C D )
1184e69d4d44SBarry Smith      the preconditioner is
1185e69d4d44SBarry Smith               (I   -B inv(A)) ( inv(A)   0    ) (I         0  )
1186e69d4d44SBarry Smith               (0    I       ) (   0    inv(S) ) (-C inv(A) I  )
1187edf189efSBarry Smith      where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1188e69d4d44SBarry Smith      inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is
1189edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1190e69d4d44SBarry Smith      D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1191e69d4d44SBarry Smith      option.
1192e69d4d44SBarry Smith 
1193edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1194edf189efSBarry Smith      is used automatically for a second block.
1195edf189efSBarry Smith 
1196a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
11970971522cSBarry Smith 
1198a541d17aSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners
1199e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
12000971522cSBarry Smith M*/
12010971522cSBarry Smith 
12020971522cSBarry Smith EXTERN_C_BEGIN
12030971522cSBarry Smith #undef __FUNCT__
12040971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
1205dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
12060971522cSBarry Smith {
12070971522cSBarry Smith   PetscErrorCode ierr;
12080971522cSBarry Smith   PC_FieldSplit  *jac;
12090971522cSBarry Smith 
12100971522cSBarry Smith   PetscFunctionBegin;
121138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
12120971522cSBarry Smith   jac->bs        = -1;
12130971522cSBarry Smith   jac->nsplits   = 0;
12143e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1215e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
121651f519a2SBarry Smith 
12170971522cSBarry Smith   pc->data     = (void*)jac;
12180971522cSBarry Smith 
12190971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1220421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
12210971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
12220971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
12230971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
12240971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
12250971522cSBarry Smith   pc->ops->applyrichardson   = 0;
12260971522cSBarry Smith 
122769a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
122869a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
12290971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
12300971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1231704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1232704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
123379416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
123479416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
123551f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
123651f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
12370971522cSBarry Smith   PetscFunctionReturn(0);
12380971522cSBarry Smith }
12390971522cSBarry Smith EXTERN_C_END
12400971522cSBarry Smith 
12410971522cSBarry Smith 
1242a541d17aSBarry Smith /*MC
1243a541d17aSBarry Smith    Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an
1244a541d17aSBarry Smith           overview of these methods.
1245a541d17aSBarry Smith 
1246a541d17aSBarry Smith       Consider the solution to ( A B ) (x_1)  =  (b_1)
1247a541d17aSBarry Smith                                ( C D ) (x_2)     (b_2)
1248a541d17aSBarry Smith 
1249a541d17aSBarry Smith       Important special cases, the Stokes equation: C = B' and D = 0  (A   B) (x_1) = (b_1)
1250a541d17aSBarry Smith                                                                        B'  0) (x_2)   (b_2)
1251a541d17aSBarry Smith 
1252a541d17aSBarry Smith       One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners
1253a541d17aSBarry Smith       for this block system.
1254a541d17aSBarry Smith 
1255a541d17aSBarry Smith       Consider an additional matrix (Ap  Bp)
1256a541d17aSBarry Smith                                     (Cp  Dp) where some or all of the entries may be the same as
1257a541d17aSBarry Smith       in the original matrix (for example Ap == A).
1258a541d17aSBarry Smith 
1259a541d17aSBarry Smith       In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the
1260a541d17aSBarry Smith       approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...)
1261a541d17aSBarry Smith 
1262a541d17aSBarry Smith       Block Jacobi:   x_1 = A^ b_1
1263a541d17aSBarry Smith                       x_2 = D^ b_2
1264a541d17aSBarry Smith 
1265a541d17aSBarry Smith       Lower block Gauss-Seidel:   x_1 = A^ b_1
1266a541d17aSBarry Smith                             x_2 = D^ (b_2 - C x_1)       variant x_2 = D^ (b_2 - Cp x_1)
1267a541d17aSBarry Smith 
1268a541d17aSBarry Smith       Symmetric Gauss-Seidel:  x_1 = x_1 + A^(b_1 - A x_1 - B x_2)    variant  x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2)
1269a541d17aSBarry Smith           Interestingly this form is not actually a symmetric matrix, the symmetric version is
1270a541d17aSBarry Smith                               x_1 = A^(b_1 - B x_2)      variant x_1 = A^(b_1 - Bp x_2)
1271a541d17aSBarry Smith 
12720bc0a719SBarry Smith    Level: intermediate
12730bc0a719SBarry Smith 
1274a541d17aSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT
1275a541d17aSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1276a541d17aSBarry Smith M*/
1277