xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 7087cfbefd1a42b179f217f9994fb6cb0d0c1824)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
36356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
40971522cSBarry Smith 
5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0};
7c5d2311dSJed Brown 
8c5d2311dSJed Brown typedef enum {
9c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG,
10c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER,
11c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER,
12c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL
13c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType;
14084e4875SJed Brown 
150971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
160971522cSBarry Smith struct _PC_FieldSplitLink {
1769a612a9SBarry Smith   KSP               ksp;
180971522cSBarry Smith   Vec               x,y;
19db4c96c1SJed Brown   char              *splitname;
200971522cSBarry Smith   PetscInt          nfields;
210971522cSBarry Smith   PetscInt          *fields;
221b9fc7fcSBarry Smith   VecScatter        sctx;
234aa3045dSJed Brown   IS                is;
2451f519a2SBarry Smith   PC_FieldSplitLink next,previous;
250971522cSBarry Smith };
260971522cSBarry Smith 
270971522cSBarry Smith typedef struct {
2868dd23aaSBarry Smith   PCCompositeType   type;
29ace3abfcSBarry Smith   PetscBool         defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
30ace3abfcSBarry Smith   PetscBool         splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
31ace3abfcSBarry Smith   PetscBool         realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
3230ad9308SMatthew Knepley   PetscInt          bs;           /* Block size for IS and Mat structures */
3330ad9308SMatthew Knepley   PetscInt          nsplits;      /* Number of field divisions defined */
3479416396SBarry Smith   Vec               *x,*y,w1,w2;
35519d70e2SJed Brown   Mat               *mat;         /* The diagonal block for each split */
36519d70e2SJed Brown   Mat               *pmat;        /* The preconditioning diagonal block for each split */
3730ad9308SMatthew Knepley   Mat               *Afield;      /* The rows of the matrix associated with each split */
38ace3abfcSBarry Smith   PetscBool         issetup;
3930ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4030ad9308SMatthew Knepley   Mat               B;            /* The (0,1) block */
4130ad9308SMatthew Knepley   Mat               C;            /* The (1,0) block */
4230ad9308SMatthew Knepley   Mat               schur;        /* The Schur complement S = D - C A^{-1} B */
43084e4875SJed Brown   Mat               schur_user;   /* User-provided preconditioning matrix for the Schur complement */
44084e4875SJed Brown   PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
45c5d2311dSJed Brown   PCFieldSplitSchurFactorizationType schurfactorization;
4630ad9308SMatthew Knepley   KSP               kspschur;     /* The solver for S */
4797bbdb24SBarry Smith   PC_FieldSplitLink head;
480971522cSBarry Smith } PC_FieldSplit;
490971522cSBarry Smith 
5016913363SBarry Smith /*
5116913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5216913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5316913363SBarry Smith    PC you could change this.
5416913363SBarry Smith */
55084e4875SJed Brown 
56e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
57084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
58084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
59084e4875SJed Brown {
60084e4875SJed Brown   switch (jac->schurpre) {
61084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
62084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
63084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
64084e4875SJed Brown     default:
65084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
66084e4875SJed Brown   }
67084e4875SJed Brown }
68084e4875SJed Brown 
69084e4875SJed Brown 
700971522cSBarry Smith #undef __FUNCT__
710971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
720971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
730971522cSBarry Smith {
740971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
750971522cSBarry Smith   PetscErrorCode    ierr;
76ace3abfcSBarry Smith   PetscBool         iascii;
770971522cSBarry Smith   PetscInt          i,j;
785a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
790971522cSBarry Smith 
800971522cSBarry Smith   PetscFunctionBegin;
812692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
820971522cSBarry Smith   if (iascii) {
8351f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8469a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
850971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
860971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
871ab39975SBarry Smith       if (ilink->fields) {
880971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
8979416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
905a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9179416396SBarry Smith 	  if (j > 0) {
9279416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9379416396SBarry Smith 	  }
945a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
950971522cSBarry Smith 	}
960971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9779416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
981ab39975SBarry Smith       } else {
991ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1001ab39975SBarry Smith       }
1015a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1025a9f2f41SSatish Balay       ilink = ilink->next;
1030971522cSBarry Smith     }
1040971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1050971522cSBarry Smith   } else {
10665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1070971522cSBarry Smith   }
1080971522cSBarry Smith   PetscFunctionReturn(0);
1090971522cSBarry Smith }
1100971522cSBarry Smith 
1110971522cSBarry Smith #undef __FUNCT__
1123b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1133b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1143b224e63SBarry Smith {
1153b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1163b224e63SBarry Smith   PetscErrorCode    ierr;
117ace3abfcSBarry Smith   PetscBool         iascii;
1183b224e63SBarry Smith   PetscInt          i,j;
1193b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1203b224e63SBarry Smith   KSP               ksp;
1213b224e63SBarry Smith 
1223b224e63SBarry Smith   PetscFunctionBegin;
1232692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1243b224e63SBarry Smith   if (iascii) {
125c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1263b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1273b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1283b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1293b224e63SBarry Smith       if (ilink->fields) {
1303b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1313b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1323b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1333b224e63SBarry Smith 	  if (j > 0) {
1343b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1353b224e63SBarry Smith 	  }
1363b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1373b224e63SBarry Smith 	}
1383b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1393b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1403b224e63SBarry Smith       } else {
1413b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1423b224e63SBarry Smith       }
1433b224e63SBarry Smith       ilink = ilink->next;
1443b224e63SBarry Smith     }
1453b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
1463b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14712cae6f2SJed Brown     if (jac->schur) {
14812cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1493b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15012cae6f2SJed Brown     } else {
15112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15212cae6f2SJed Brown     }
1533b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1543b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
1553b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15612cae6f2SJed Brown     if (jac->kspschur) {
1573b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
15812cae6f2SJed Brown     } else {
15912cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16012cae6f2SJed Brown     }
1613b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1623b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1633b224e63SBarry Smith   } else {
16465e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1653b224e63SBarry Smith   }
1663b224e63SBarry Smith   PetscFunctionReturn(0);
1673b224e63SBarry Smith }
1683b224e63SBarry Smith 
1693b224e63SBarry Smith #undef __FUNCT__
1706c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1716c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1726c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1736c924f48SJed Brown {
1746c924f48SJed Brown   PetscErrorCode ierr;
1756c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1766c924f48SJed Brown   PetscInt       i,nfields,*ifields;
177ace3abfcSBarry Smith   PetscBool      flg;
1786c924f48SJed Brown   char           optionname[128],splitname[8];
1796c924f48SJed Brown 
1806c924f48SJed Brown   PetscFunctionBegin;
1816c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1826c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1836c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1846c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1856c924f48SJed Brown     nfields = jac->bs;
1866c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1876c924f48SJed Brown     if (!flg) break;
1886c924f48SJed Brown     if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1896c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1906c924f48SJed Brown   }
1916c924f48SJed Brown   if (i > 0) {
1926c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1936c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1946c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1956c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1966c924f48SJed Brown   }
1976c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
1986c924f48SJed Brown   PetscFunctionReturn(0);
1996c924f48SJed Brown }
2006c924f48SJed Brown 
2016c924f48SJed Brown #undef __FUNCT__
20269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2040971522cSBarry Smith {
2050971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2060971522cSBarry Smith   PetscErrorCode    ierr;
2075a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
208ace3abfcSBarry Smith   PetscBool         flg = PETSC_FALSE;
2096c924f48SJed Brown   PetscInt          i;
2100971522cSBarry Smith 
2110971522cSBarry Smith   PetscFunctionBegin;
212d32f9abdSBarry Smith   if (!ilink) {
2138b8307b2SJed Brown     if (pc->dm) {
2148b8307b2SJed Brown       PetscBool dmcomposite;
2158b8307b2SJed Brown       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
2168b8307b2SJed Brown       if (dmcomposite) {
2178b8307b2SJed Brown         PetscInt nDM;
2188b8307b2SJed Brown         IS       *fields;
2198b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2208b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
2218b8307b2SJed Brown         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
2228b8307b2SJed Brown         for (i=0; i<nDM; i++) {
2238b8307b2SJed Brown           char splitname[8];
2248b8307b2SJed Brown           ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2258b8307b2SJed Brown           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
2268b8307b2SJed Brown           ierr = ISDestroy(fields[i]);CHKERRQ(ierr);
2278b8307b2SJed Brown         }
2288b8307b2SJed Brown         ierr = PetscFree(fields);CHKERRQ(ierr);
2298b8307b2SJed Brown       }
23066ffff09SJed Brown     } else {
231521d7252SBarry Smith       if (jac->bs <= 0) {
232704ba839SBarry Smith         if (pc->pmat) {
233521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
234704ba839SBarry Smith         } else {
235704ba839SBarry Smith           jac->bs = 1;
236704ba839SBarry Smith         }
237521d7252SBarry Smith       }
238d32f9abdSBarry Smith 
239acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
240d32f9abdSBarry Smith       if (!flg) {
241d32f9abdSBarry Smith         /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
242d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2436c924f48SJed Brown         ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2446c924f48SJed Brown         if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
245d32f9abdSBarry Smith       }
2466c924f48SJed Brown       if (flg || !jac->splitdefined) {
247d32f9abdSBarry Smith         ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
248db4c96c1SJed Brown         for (i=0; i<jac->bs; i++) {
2496c924f48SJed Brown           char splitname[8];
2506c924f48SJed Brown           ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
251db4c96c1SJed Brown           ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
25279416396SBarry Smith         }
25397bbdb24SBarry Smith         jac->defaultsplit = PETSC_TRUE;
254521d7252SBarry Smith       }
25566ffff09SJed Brown     }
256edf189efSBarry Smith   } else if (jac->nsplits == 1) {
257edf189efSBarry Smith     if (ilink->is) {
258edf189efSBarry Smith       IS       is2;
259edf189efSBarry Smith       PetscInt nmin,nmax;
260edf189efSBarry Smith 
261edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
262edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
263db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
264edf189efSBarry Smith       ierr = ISDestroy(is2);CHKERRQ(ierr);
265e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
266edf189efSBarry Smith   }
267e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
26869a612a9SBarry Smith   PetscFunctionReturn(0);
26969a612a9SBarry Smith }
27069a612a9SBarry Smith 
27169a612a9SBarry Smith 
27269a612a9SBarry Smith #undef __FUNCT__
27369a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
27469a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
27569a612a9SBarry Smith {
27669a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
27769a612a9SBarry Smith   PetscErrorCode    ierr;
2785a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
279cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
28069a612a9SBarry Smith   MatStructure      flag = pc->flag;
281ace3abfcSBarry Smith   PetscBool         sorted;
28269a612a9SBarry Smith 
28369a612a9SBarry Smith   PetscFunctionBegin;
28469a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
28597bbdb24SBarry Smith   nsplit = jac->nsplits;
2865a9f2f41SSatish Balay   ilink  = jac->head;
28797bbdb24SBarry Smith 
28897bbdb24SBarry Smith   /* get the matrices for each split */
289704ba839SBarry Smith   if (!jac->issetup) {
2901b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
29197bbdb24SBarry Smith 
292704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
293704ba839SBarry Smith 
294704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
29551f519a2SBarry Smith     bs     = jac->bs;
29697bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
297cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
2981b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
2991b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3001b9fc7fcSBarry Smith       if (jac->defaultsplit) {
301704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
302704ba839SBarry Smith       } else if (!ilink->is) {
303ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3045a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3055a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3061b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3071b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3081b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
30997bbdb24SBarry Smith             }
31097bbdb24SBarry Smith           }
311d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
312ccb205f8SBarry Smith         } else {
313704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
314ccb205f8SBarry Smith         }
3153e197d65SBarry Smith       }
316704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
317e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
318704ba839SBarry Smith       ilink = ilink->next;
3191b9fc7fcSBarry Smith     }
3201b9fc7fcSBarry Smith   }
3211b9fc7fcSBarry Smith 
322704ba839SBarry Smith   ilink  = jac->head;
32397bbdb24SBarry Smith   if (!jac->pmat) {
324cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
325cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3264aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
327704ba839SBarry Smith       ilink = ilink->next;
328cf502942SBarry Smith     }
32997bbdb24SBarry Smith   } else {
330cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3314aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
332704ba839SBarry Smith       ilink = ilink->next;
333cf502942SBarry Smith     }
33497bbdb24SBarry Smith   }
335519d70e2SJed Brown   if (jac->realdiagonal) {
336519d70e2SJed Brown     ilink = jac->head;
337519d70e2SJed Brown     if (!jac->mat) {
338519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
339519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
340519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
341519d70e2SJed Brown         ilink = ilink->next;
342519d70e2SJed Brown       }
343519d70e2SJed Brown     } else {
344519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
345519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
346519d70e2SJed Brown         ilink = ilink->next;
347519d70e2SJed Brown       }
348519d70e2SJed Brown     }
349519d70e2SJed Brown   } else {
350519d70e2SJed Brown     jac->mat = jac->pmat;
351519d70e2SJed Brown   }
35297bbdb24SBarry Smith 
3536c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
35468dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
35568dd23aaSBarry Smith     ilink  = jac->head;
35668dd23aaSBarry Smith     if (!jac->Afield) {
35768dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
35868dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3594aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
36068dd23aaSBarry Smith         ilink = ilink->next;
36168dd23aaSBarry Smith       }
36268dd23aaSBarry Smith     } else {
36368dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3644aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
36568dd23aaSBarry Smith         ilink = ilink->next;
36668dd23aaSBarry Smith       }
36768dd23aaSBarry Smith     }
36868dd23aaSBarry Smith   }
36968dd23aaSBarry Smith 
3703b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
3713b224e63SBarry Smith     IS       ccis;
3724aa3045dSJed Brown     PetscInt rstart,rend;
373e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
37468dd23aaSBarry Smith 
375e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
376e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
377e6cab6aaSJed Brown 
3783b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
3793b224e63SBarry Smith     if (jac->schur) {
3803b224e63SBarry Smith       ilink = jac->head;
3814aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3824aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
3833b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3843b224e63SBarry Smith       ilink = ilink->next;
3854aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3864aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
3873b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
388519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
389084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
3903b224e63SBarry Smith 
3913b224e63SBarry Smith      } else {
3921cee3971SBarry Smith       KSP ksp;
3936c924f48SJed Brown       char schurprefix[256];
3943b224e63SBarry Smith 
3953b224e63SBarry Smith       /* extract the B and C matrices */
3963b224e63SBarry Smith       ilink = jac->head;
3974aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3984aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
3993b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
4003b224e63SBarry Smith       ilink = ilink->next;
4014aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4024aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
4033b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
404084e4875SJed Brown       /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */
405519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
4061cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
407e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
4083b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4093b224e63SBarry Smith 
4103b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4119005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4121cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
413084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
414084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
415084e4875SJed Brown         PC pc;
416cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
417084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
418084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
419e69d4d44SBarry Smith       }
4206c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4216c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4223b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4233b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4243b224e63SBarry Smith 
4253b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4263b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4273b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4283b224e63SBarry Smith       ilink = jac->head;
4293b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4303b224e63SBarry Smith       ilink = ilink->next;
4313b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4323b224e63SBarry Smith     }
4333b224e63SBarry Smith   } else {
43497bbdb24SBarry Smith     /* set up the individual PCs */
43597bbdb24SBarry Smith     i    = 0;
4365a9f2f41SSatish Balay     ilink = jac->head;
4375a9f2f41SSatish Balay     while (ilink) {
438519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
4393b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4405a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
4415a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
44297bbdb24SBarry Smith       i++;
4435a9f2f41SSatish Balay       ilink = ilink->next;
4440971522cSBarry Smith     }
44597bbdb24SBarry Smith 
44697bbdb24SBarry Smith     /* create work vectors for each split */
4471b9fc7fcSBarry Smith     if (!jac->x) {
44897bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
4495a9f2f41SSatish Balay       ilink = jac->head;
45097bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
451906ed7ccSBarry Smith         Vec *vl,*vr;
452906ed7ccSBarry Smith 
453906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
454906ed7ccSBarry Smith         ilink->x  = *vr;
455906ed7ccSBarry Smith         ilink->y  = *vl;
456906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
457906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
4585a9f2f41SSatish Balay         jac->x[i] = ilink->x;
4595a9f2f41SSatish Balay         jac->y[i] = ilink->y;
4605a9f2f41SSatish Balay         ilink     = ilink->next;
46197bbdb24SBarry Smith       }
4623b224e63SBarry Smith     }
4633b224e63SBarry Smith   }
4643b224e63SBarry Smith 
4653b224e63SBarry Smith 
466d0f46423SBarry Smith   if (!jac->head->sctx) {
4673b224e63SBarry Smith     Vec xtmp;
4683b224e63SBarry Smith 
46979416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
4701b9fc7fcSBarry Smith 
4715a9f2f41SSatish Balay     ilink = jac->head;
4721b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
4731b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
474704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
4755a9f2f41SSatish Balay       ilink = ilink->next;
4761b9fc7fcSBarry Smith     }
4771b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
4781b9fc7fcSBarry Smith   }
4790971522cSBarry Smith   PetscFunctionReturn(0);
4800971522cSBarry Smith }
4810971522cSBarry Smith 
4825a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
483ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
484ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
4855a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
486ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
487ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
48879416396SBarry Smith 
4890971522cSBarry Smith #undef __FUNCT__
4903b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
4913b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
4923b224e63SBarry Smith {
4933b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4943b224e63SBarry Smith   PetscErrorCode    ierr;
4953b224e63SBarry Smith   KSP               ksp;
4963b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
4973b224e63SBarry Smith 
4983b224e63SBarry Smith   PetscFunctionBegin;
4993b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5003b224e63SBarry Smith 
501c5d2311dSJed Brown   switch (jac->schurfactorization) {
502c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
503c5d2311dSJed Brown       /* [A 0; 0 -S], positive definite, suitable for MINRES */
504c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
505c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
506c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
507c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
508c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
509c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
510c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
511c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
512c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
513c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
514c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
515c5d2311dSJed Brown       break;
516c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
517c5d2311dSJed Brown       /* [A 0; C S], suitable for left preconditioning */
518c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
519c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
520c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
521c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
522c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
523c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
524c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
525c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
526c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
527c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
528c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
529c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530c5d2311dSJed Brown       break;
531c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
532c5d2311dSJed Brown       /* [A B; 0 S], suitable for right preconditioning */
533c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
534c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
535c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
536c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
537c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
538c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
539c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
540c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
541c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
542c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
543c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
544c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
545c5d2311dSJed Brown       break;
546c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
547c5d2311dSJed Brown       /* [1 0; CA^{-1} 1] [A 0; 0 S] [1 A^{-1}B; 0 1], an exact solve if applied exactly, needs one extra solve with A */
5483b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5493b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5503b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5513b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
5523b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
5533b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5543b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5553b224e63SBarry Smith 
5563b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
5573b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5583b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5593b224e63SBarry Smith 
5603b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
5613b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
5623b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5633b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5643b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
565c5d2311dSJed Brown   }
5663b224e63SBarry Smith   PetscFunctionReturn(0);
5673b224e63SBarry Smith }
5683b224e63SBarry Smith 
5693b224e63SBarry Smith #undef __FUNCT__
5700971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
5710971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
5720971522cSBarry Smith {
5730971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5740971522cSBarry Smith   PetscErrorCode    ierr;
5755a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
576af93645bSJed Brown   PetscInt          cnt;
5770971522cSBarry Smith 
5780971522cSBarry Smith   PetscFunctionBegin;
57951f519a2SBarry Smith   CHKMEMQ;
58051f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
58151f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
58251f519a2SBarry Smith 
58379416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
5841b9fc7fcSBarry Smith     if (jac->defaultsplit) {
5850971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
5865a9f2f41SSatish Balay       while (ilink) {
5875a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
5885a9f2f41SSatish Balay         ilink = ilink->next;
5890971522cSBarry Smith       }
5900971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
5911b9fc7fcSBarry Smith     } else {
592efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
5935a9f2f41SSatish Balay       while (ilink) {
5945a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
5955a9f2f41SSatish Balay         ilink = ilink->next;
5961b9fc7fcSBarry Smith       }
5971b9fc7fcSBarry Smith     }
59816913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
59979416396SBarry Smith     if (!jac->w1) {
60079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
60179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
60279416396SBarry Smith     }
603efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6045a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6053e197d65SBarry Smith     cnt = 1;
6065a9f2f41SSatish Balay     while (ilink->next) {
6075a9f2f41SSatish Balay       ilink = ilink->next;
6083e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6093e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6103e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6113e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6123e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6133e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6143e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6153e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6163e197d65SBarry Smith     }
61751f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
61811755939SBarry Smith       cnt -= 2;
61951f519a2SBarry Smith       while (ilink->previous) {
62051f519a2SBarry Smith         ilink = ilink->previous;
62111755939SBarry Smith         /* compute the residual only over the part of the vector needed */
62211755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
62311755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
62411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62611755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
62711755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
62811755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
62951f519a2SBarry Smith       }
63011755939SBarry Smith     }
63165e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
63251f519a2SBarry Smith   CHKMEMQ;
6330971522cSBarry Smith   PetscFunctionReturn(0);
6340971522cSBarry Smith }
6350971522cSBarry Smith 
636421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
637ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
638ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
639421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
640ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
641ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
642421e10b8SBarry Smith 
643421e10b8SBarry Smith #undef __FUNCT__
6448c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
645421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
646421e10b8SBarry Smith {
647421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
648421e10b8SBarry Smith   PetscErrorCode    ierr;
649421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
650421e10b8SBarry Smith 
651421e10b8SBarry Smith   PetscFunctionBegin;
652421e10b8SBarry Smith   CHKMEMQ;
653421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
654421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
655421e10b8SBarry Smith 
656421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
657421e10b8SBarry Smith     if (jac->defaultsplit) {
658421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
659421e10b8SBarry Smith       while (ilink) {
660421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
661421e10b8SBarry Smith 	ilink = ilink->next;
662421e10b8SBarry Smith       }
663421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
664421e10b8SBarry Smith     } else {
665421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
666421e10b8SBarry Smith       while (ilink) {
667421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
668421e10b8SBarry Smith 	ilink = ilink->next;
669421e10b8SBarry Smith       }
670421e10b8SBarry Smith     }
671421e10b8SBarry Smith   } else {
672421e10b8SBarry Smith     if (!jac->w1) {
673421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
674421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
675421e10b8SBarry Smith     }
676421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
677421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
678421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
679421e10b8SBarry Smith       while (ilink->next) {
680421e10b8SBarry Smith         ilink = ilink->next;
6819989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
682421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
683421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
684421e10b8SBarry Smith       }
685421e10b8SBarry Smith       while (ilink->previous) {
686421e10b8SBarry Smith         ilink = ilink->previous;
6879989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
688421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
689421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
690421e10b8SBarry Smith       }
691421e10b8SBarry Smith     } else {
692421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
693421e10b8SBarry Smith 	ilink = ilink->next;
694421e10b8SBarry Smith       }
695421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
696421e10b8SBarry Smith       while (ilink->previous) {
697421e10b8SBarry Smith 	ilink = ilink->previous;
6989989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
699421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
700421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
701421e10b8SBarry Smith       }
702421e10b8SBarry Smith     }
703421e10b8SBarry Smith   }
704421e10b8SBarry Smith   CHKMEMQ;
705421e10b8SBarry Smith   PetscFunctionReturn(0);
706421e10b8SBarry Smith }
707421e10b8SBarry Smith 
7080971522cSBarry Smith #undef __FUNCT__
7090971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
7100971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
7110971522cSBarry Smith {
7120971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7130971522cSBarry Smith   PetscErrorCode    ierr;
7145a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7150971522cSBarry Smith 
7160971522cSBarry Smith   PetscFunctionBegin;
7175a9f2f41SSatish Balay   while (ilink) {
7185a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
7195a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
7205a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
7215a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
722704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
7235a9f2f41SSatish Balay     next = ilink->next;
7247e8c30b6SJed Brown     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
725704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
726704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
7275a9f2f41SSatish Balay     ilink = next;
7280971522cSBarry Smith   }
72905b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
730519d70e2SJed Brown   if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);}
73197bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
73268dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
73379416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
73479416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
7353b224e63SBarry Smith   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
736084e4875SJed Brown   if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
737d0f46423SBarry Smith   if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);}
7383b224e63SBarry Smith   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
7393b224e63SBarry Smith   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
7400971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
7410971522cSBarry Smith   PetscFunctionReturn(0);
7420971522cSBarry Smith }
7430971522cSBarry Smith 
7440971522cSBarry Smith #undef __FUNCT__
7450971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
7460971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
7470971522cSBarry Smith {
7481b9fc7fcSBarry Smith   PetscErrorCode  ierr;
7496c924f48SJed Brown   PetscInt        bs;
750ace3abfcSBarry Smith   PetscBool       flg;
7519dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
7523b224e63SBarry Smith   PCCompositeType ctype;
7531b9fc7fcSBarry Smith 
7540971522cSBarry Smith   PetscFunctionBegin;
7551b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
756acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
75751f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
75851f519a2SBarry Smith   if (flg) {
75951f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
76051f519a2SBarry Smith   }
761704ba839SBarry Smith 
7623b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
7633b224e63SBarry Smith   if (flg) {
7643b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
7653b224e63SBarry Smith   }
766d32f9abdSBarry Smith 
767c30613efSMatthew Knepley   /* Only setup fields once */
768c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
769d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
770d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
7716c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
7726c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
773d32f9abdSBarry Smith   }
774c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
775c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
776084e4875SJed 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);
777c5d2311dSJed Brown   }
7781b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7790971522cSBarry Smith   PetscFunctionReturn(0);
7800971522cSBarry Smith }
7810971522cSBarry Smith 
7820971522cSBarry Smith /*------------------------------------------------------------------------------------*/
7830971522cSBarry Smith 
7840971522cSBarry Smith EXTERN_C_BEGIN
7850971522cSBarry Smith #undef __FUNCT__
7860971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
787*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
7880971522cSBarry Smith {
78997bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7900971522cSBarry Smith   PetscErrorCode    ierr;
7915a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
79269a612a9SBarry Smith   char              prefix[128];
79351f519a2SBarry Smith   PetscInt          i;
7940971522cSBarry Smith 
7950971522cSBarry Smith   PetscFunctionBegin;
7966c924f48SJed Brown   if (jac->splitdefined) {
7976c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
7986c924f48SJed Brown     PetscFunctionReturn(0);
7996c924f48SJed Brown   }
80051f519a2SBarry Smith   for (i=0; i<n; i++) {
801e32f2f54SBarry 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);
802e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
80351f519a2SBarry Smith   }
804704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
805db4c96c1SJed Brown   ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
806704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
8075a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
8085a9f2f41SSatish Balay   ilink->nfields = n;
8095a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
8107adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
8111cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
8125a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
8139005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
81469a612a9SBarry Smith 
8156c924f48SJed Brown   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr);
8165a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
8170971522cSBarry Smith 
8180971522cSBarry Smith   if (!next) {
8195a9f2f41SSatish Balay     jac->head       = ilink;
82051f519a2SBarry Smith     ilink->previous = PETSC_NULL;
8210971522cSBarry Smith   } else {
8220971522cSBarry Smith     while (next->next) {
8230971522cSBarry Smith       next = next->next;
8240971522cSBarry Smith     }
8255a9f2f41SSatish Balay     next->next      = ilink;
82651f519a2SBarry Smith     ilink->previous = next;
8270971522cSBarry Smith   }
8280971522cSBarry Smith   jac->nsplits++;
8290971522cSBarry Smith   PetscFunctionReturn(0);
8300971522cSBarry Smith }
8310971522cSBarry Smith EXTERN_C_END
8320971522cSBarry Smith 
833e69d4d44SBarry Smith EXTERN_C_BEGIN
834e69d4d44SBarry Smith #undef __FUNCT__
835e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
836*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
837e69d4d44SBarry Smith {
838e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
839e69d4d44SBarry Smith   PetscErrorCode ierr;
840e69d4d44SBarry Smith 
841e69d4d44SBarry Smith   PetscFunctionBegin;
842519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
843e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
844e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
845084e4875SJed Brown   *n = jac->nsplits;
846e69d4d44SBarry Smith   PetscFunctionReturn(0);
847e69d4d44SBarry Smith }
848e69d4d44SBarry Smith EXTERN_C_END
8490971522cSBarry Smith 
8500971522cSBarry Smith EXTERN_C_BEGIN
8510971522cSBarry Smith #undef __FUNCT__
85269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
853*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
8540971522cSBarry Smith {
8550971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8560971522cSBarry Smith   PetscErrorCode    ierr;
8570971522cSBarry Smith   PetscInt          cnt = 0;
8585a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
8590971522cSBarry Smith 
8600971522cSBarry Smith   PetscFunctionBegin;
86169a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
8625a9f2f41SSatish Balay   while (ilink) {
8635a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
8645a9f2f41SSatish Balay     ilink = ilink->next;
8650971522cSBarry Smith   }
866e32f2f54SBarry 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);
8670971522cSBarry Smith   *n = jac->nsplits;
8680971522cSBarry Smith   PetscFunctionReturn(0);
8690971522cSBarry Smith }
8700971522cSBarry Smith EXTERN_C_END
8710971522cSBarry Smith 
872704ba839SBarry Smith EXTERN_C_BEGIN
873704ba839SBarry Smith #undef __FUNCT__
874704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
875*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
876704ba839SBarry Smith {
877704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
878704ba839SBarry Smith   PetscErrorCode    ierr;
879704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
880704ba839SBarry Smith   char              prefix[128];
881704ba839SBarry Smith 
882704ba839SBarry Smith   PetscFunctionBegin;
8836c924f48SJed Brown   if (jac->splitdefined) {
8846c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8856c924f48SJed Brown     PetscFunctionReturn(0);
8866c924f48SJed Brown   }
88716913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
888db4c96c1SJed Brown   ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
8891ab39975SBarry Smith   ilink->is      = is;
8901ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
891704ba839SBarry Smith   ilink->next    = PETSC_NULL;
892704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
8931cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
894704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
8959005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
896704ba839SBarry Smith 
8976c924f48SJed Brown   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr);
898704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
899704ba839SBarry Smith 
900704ba839SBarry Smith   if (!next) {
901704ba839SBarry Smith     jac->head       = ilink;
902704ba839SBarry Smith     ilink->previous = PETSC_NULL;
903704ba839SBarry Smith   } else {
904704ba839SBarry Smith     while (next->next) {
905704ba839SBarry Smith       next = next->next;
906704ba839SBarry Smith     }
907704ba839SBarry Smith     next->next      = ilink;
908704ba839SBarry Smith     ilink->previous = next;
909704ba839SBarry Smith   }
910704ba839SBarry Smith   jac->nsplits++;
911704ba839SBarry Smith 
912704ba839SBarry Smith   PetscFunctionReturn(0);
913704ba839SBarry Smith }
914704ba839SBarry Smith EXTERN_C_END
915704ba839SBarry Smith 
9160971522cSBarry Smith #undef __FUNCT__
9170971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
9180971522cSBarry Smith /*@
9190971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
9200971522cSBarry Smith 
921ad4df100SBarry Smith     Logically Collective on PC
9220971522cSBarry Smith 
9230971522cSBarry Smith     Input Parameters:
9240971522cSBarry Smith +   pc  - the preconditioner context
925db4c96c1SJed Brown .   splitname - name of this split
9260971522cSBarry Smith .   n - the number of fields in this split
927db4c96c1SJed Brown -   fields - the fields in this split
9280971522cSBarry Smith 
9290971522cSBarry Smith     Level: intermediate
9300971522cSBarry Smith 
931d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
932d32f9abdSBarry Smith 
933d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
934d32f9abdSBarry 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
935d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
936d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
937d32f9abdSBarry Smith 
938db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
939db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
940db4c96c1SJed Brown 
941d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
9420971522cSBarry Smith 
9430971522cSBarry Smith @*/
944*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
9450971522cSBarry Smith {
9464ac538c5SBarry Smith   PetscErrorCode ierr;
9470971522cSBarry Smith 
9480971522cSBarry Smith   PetscFunctionBegin;
9490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
950db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
951db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
952db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
9534ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
9540971522cSBarry Smith   PetscFunctionReturn(0);
9550971522cSBarry Smith }
9560971522cSBarry Smith 
9570971522cSBarry Smith #undef __FUNCT__
958704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
959704ba839SBarry Smith /*@
960704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
961704ba839SBarry Smith 
962ad4df100SBarry Smith     Logically Collective on PC
963704ba839SBarry Smith 
964704ba839SBarry Smith     Input Parameters:
965704ba839SBarry Smith +   pc  - the preconditioner context
966db4c96c1SJed Brown .   splitname - name of this split
967db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
968704ba839SBarry Smith 
969d32f9abdSBarry Smith 
970a6ffb8dbSJed Brown     Notes:
971a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
972a6ffb8dbSJed Brown 
973db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
974db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
975d32f9abdSBarry Smith 
976704ba839SBarry Smith     Level: intermediate
977704ba839SBarry Smith 
978704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
979704ba839SBarry Smith 
980704ba839SBarry Smith @*/
981*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
982704ba839SBarry Smith {
9834ac538c5SBarry Smith   PetscErrorCode ierr;
984704ba839SBarry Smith 
985704ba839SBarry Smith   PetscFunctionBegin;
9860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
987db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
988db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
9894ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
990704ba839SBarry Smith   PetscFunctionReturn(0);
991704ba839SBarry Smith }
992704ba839SBarry Smith 
993704ba839SBarry Smith #undef __FUNCT__
99451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
99551f519a2SBarry Smith /*@
99651f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
99751f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
99851f519a2SBarry Smith 
999ad4df100SBarry Smith     Logically Collective on PC
100051f519a2SBarry Smith 
100151f519a2SBarry Smith     Input Parameters:
100251f519a2SBarry Smith +   pc  - the preconditioner context
100351f519a2SBarry Smith -   bs - the block size
100451f519a2SBarry Smith 
100551f519a2SBarry Smith     Level: intermediate
100651f519a2SBarry Smith 
100751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
100851f519a2SBarry Smith 
100951f519a2SBarry Smith @*/
1010*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
101151f519a2SBarry Smith {
10124ac538c5SBarry Smith   PetscErrorCode ierr;
101351f519a2SBarry Smith 
101451f519a2SBarry Smith   PetscFunctionBegin;
10150700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1016c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
10174ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
101851f519a2SBarry Smith   PetscFunctionReturn(0);
101951f519a2SBarry Smith }
102051f519a2SBarry Smith 
102151f519a2SBarry Smith #undef __FUNCT__
102269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
10230971522cSBarry Smith /*@C
102469a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
10250971522cSBarry Smith 
102669a612a9SBarry Smith    Collective on KSP
10270971522cSBarry Smith 
10280971522cSBarry Smith    Input Parameter:
10290971522cSBarry Smith .  pc - the preconditioner context
10300971522cSBarry Smith 
10310971522cSBarry Smith    Output Parameters:
10320971522cSBarry Smith +  n - the number of split
103369a612a9SBarry Smith -  pc - the array of KSP contexts
10340971522cSBarry Smith 
10350971522cSBarry Smith    Note:
1036d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1037d32f9abdSBarry Smith    (not the KSP just the array that contains them).
10380971522cSBarry Smith 
103969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
10400971522cSBarry Smith 
10410971522cSBarry Smith    Level: advanced
10420971522cSBarry Smith 
10430971522cSBarry Smith .seealso: PCFIELDSPLIT
10440971522cSBarry Smith @*/
1045*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
10460971522cSBarry Smith {
10474ac538c5SBarry Smith   PetscErrorCode ierr;
10480971522cSBarry Smith 
10490971522cSBarry Smith   PetscFunctionBegin;
10500700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10510971522cSBarry Smith   PetscValidIntPointer(n,2);
10524ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
10530971522cSBarry Smith   PetscFunctionReturn(0);
10540971522cSBarry Smith }
10550971522cSBarry Smith 
1056e69d4d44SBarry Smith #undef __FUNCT__
1057e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1058e69d4d44SBarry Smith /*@
1059e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1060e69d4d44SBarry Smith       D matrix. Otherwise no preconditioner is used.
1061e69d4d44SBarry Smith 
1062e69d4d44SBarry Smith     Collective on PC
1063e69d4d44SBarry Smith 
1064e69d4d44SBarry Smith     Input Parameters:
1065e69d4d44SBarry Smith +   pc  - the preconditioner context
1066084e4875SJed Brown .   ptype - which matrix to use for preconditioning the Schur complement
1067084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1068084e4875SJed Brown 
1069084e4875SJed Brown     Notes:
1070084e4875SJed Brown     The default is to use the block on the diagonal of the preconditioning matrix.  This is D, in the (1,1) position.
1071084e4875SJed Brown     There are currently no preconditioners that work directly with the Schur complement so setting
1072084e4875SJed Brown     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
1073e69d4d44SBarry Smith 
1074e69d4d44SBarry Smith     Options Database:
1075084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1076e69d4d44SBarry Smith 
1077e69d4d44SBarry Smith     Level: intermediate
1078e69d4d44SBarry Smith 
1079084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1080e69d4d44SBarry Smith 
1081e69d4d44SBarry Smith @*/
1082*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1083e69d4d44SBarry Smith {
10844ac538c5SBarry Smith   PetscErrorCode ierr;
1085e69d4d44SBarry Smith 
1086e69d4d44SBarry Smith   PetscFunctionBegin;
10870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10884ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1089e69d4d44SBarry Smith   PetscFunctionReturn(0);
1090e69d4d44SBarry Smith }
1091e69d4d44SBarry Smith 
1092e69d4d44SBarry Smith EXTERN_C_BEGIN
1093e69d4d44SBarry Smith #undef __FUNCT__
1094e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1095*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1096e69d4d44SBarry Smith {
1097e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1098084e4875SJed Brown   PetscErrorCode  ierr;
1099e69d4d44SBarry Smith 
1100e69d4d44SBarry Smith   PetscFunctionBegin;
1101084e4875SJed Brown   jac->schurpre = ptype;
1102084e4875SJed Brown   if (pre) {
1103084e4875SJed Brown     if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
1104084e4875SJed Brown     jac->schur_user = pre;
1105084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1106084e4875SJed Brown   }
1107e69d4d44SBarry Smith   PetscFunctionReturn(0);
1108e69d4d44SBarry Smith }
1109e69d4d44SBarry Smith EXTERN_C_END
1110e69d4d44SBarry Smith 
111130ad9308SMatthew Knepley #undef __FUNCT__
111230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
111330ad9308SMatthew Knepley /*@C
11148c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
111530ad9308SMatthew Knepley 
111630ad9308SMatthew Knepley    Collective on KSP
111730ad9308SMatthew Knepley 
111830ad9308SMatthew Knepley    Input Parameter:
111930ad9308SMatthew Knepley .  pc - the preconditioner context
112030ad9308SMatthew Knepley 
112130ad9308SMatthew Knepley    Output Parameters:
112230ad9308SMatthew Knepley +  A - the (0,0) block
112330ad9308SMatthew Knepley .  B - the (0,1) block
112430ad9308SMatthew Knepley .  C - the (1,0) block
112530ad9308SMatthew Knepley -  D - the (1,1) block
112630ad9308SMatthew Knepley 
112730ad9308SMatthew Knepley    Level: advanced
112830ad9308SMatthew Knepley 
112930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
113030ad9308SMatthew Knepley @*/
1131*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D)
113230ad9308SMatthew Knepley {
113330ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
113430ad9308SMatthew Knepley 
113530ad9308SMatthew Knepley   PetscFunctionBegin;
11360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1137c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
113830ad9308SMatthew Knepley   if (A) *A = jac->pmat[0];
113930ad9308SMatthew Knepley   if (B) *B = jac->B;
114030ad9308SMatthew Knepley   if (C) *C = jac->C;
114130ad9308SMatthew Knepley   if (D) *D = jac->pmat[1];
114230ad9308SMatthew Knepley   PetscFunctionReturn(0);
114330ad9308SMatthew Knepley }
114430ad9308SMatthew Knepley 
114579416396SBarry Smith EXTERN_C_BEGIN
114679416396SBarry Smith #undef __FUNCT__
114779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1148*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
114979416396SBarry Smith {
115079416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1151e69d4d44SBarry Smith   PetscErrorCode ierr;
115279416396SBarry Smith 
115379416396SBarry Smith   PetscFunctionBegin;
115479416396SBarry Smith   jac->type = type;
11553b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
11563b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
11573b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1158e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1159e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1160e69d4d44SBarry Smith 
11613b224e63SBarry Smith   } else {
11623b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
11633b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1164e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
11659e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
11663b224e63SBarry Smith   }
116779416396SBarry Smith   PetscFunctionReturn(0);
116879416396SBarry Smith }
116979416396SBarry Smith EXTERN_C_END
117079416396SBarry Smith 
117151f519a2SBarry Smith EXTERN_C_BEGIN
117251f519a2SBarry Smith #undef __FUNCT__
117351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
1174*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
117551f519a2SBarry Smith {
117651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
117751f519a2SBarry Smith 
117851f519a2SBarry Smith   PetscFunctionBegin;
117965e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
118065e19b50SBarry 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);
118151f519a2SBarry Smith   jac->bs = bs;
118251f519a2SBarry Smith   PetscFunctionReturn(0);
118351f519a2SBarry Smith }
118451f519a2SBarry Smith EXTERN_C_END
118551f519a2SBarry Smith 
118679416396SBarry Smith #undef __FUNCT__
118779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1188bc08b0f1SBarry Smith /*@
118979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
119079416396SBarry Smith 
119179416396SBarry Smith    Collective on PC
119279416396SBarry Smith 
119379416396SBarry Smith    Input Parameter:
119479416396SBarry Smith .  pc - the preconditioner context
119581540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
119679416396SBarry Smith 
119779416396SBarry Smith    Options Database Key:
1198a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
119979416396SBarry Smith 
120079416396SBarry Smith    Level: Developer
120179416396SBarry Smith 
120279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
120379416396SBarry Smith 
120479416396SBarry Smith .seealso: PCCompositeSetType()
120579416396SBarry Smith 
120679416396SBarry Smith @*/
1207*7087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
120879416396SBarry Smith {
12094ac538c5SBarry Smith   PetscErrorCode ierr;
121079416396SBarry Smith 
121179416396SBarry Smith   PetscFunctionBegin;
12120700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12134ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
121479416396SBarry Smith   PetscFunctionReturn(0);
121579416396SBarry Smith }
121679416396SBarry Smith 
12170971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
12180971522cSBarry Smith /*MC
1219a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
12200971522cSBarry Smith                   fields or groups of fields
12210971522cSBarry Smith 
12220971522cSBarry Smith 
1223edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1224edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
12250971522cSBarry Smith 
1226a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
122769a612a9SBarry Smith          and set the options directly on the resulting KSP object
12280971522cSBarry Smith 
12290971522cSBarry Smith    Level: intermediate
12300971522cSBarry Smith 
123179416396SBarry Smith    Options Database Keys:
123281540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
123381540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
123481540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
123581540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
123681540f2fSBarry Smith .    -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1237e69d4d44SBarry Smith .    -pc_fieldsplit_schur_precondition <true,false> default is true
123879416396SBarry Smith 
1239edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
12403b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
12413b224e63SBarry Smith 
12423b224e63SBarry Smith 
1243d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1244d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1245d32f9abdSBarry Smith 
1246d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1247d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1248d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1249d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1250d32f9abdSBarry Smith 
1251d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
1252d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
1253d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
1254d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1255d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
1256d32f9abdSBarry Smith 
1257e69d4d44SBarry Smith      For the Schur complement preconditioner if J = ( A B )
1258e69d4d44SBarry Smith                                                     ( C D )
1259e69d4d44SBarry Smith      the preconditioner is
1260e69d4d44SBarry Smith               (I   -B inv(A)) ( inv(A)   0    ) (I         0  )
1261e69d4d44SBarry Smith               (0    I       ) (   0    inv(S) ) (-C inv(A) I  )
1262edf189efSBarry Smith      where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1263e69d4d44SBarry Smith      inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is
1264edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1265e69d4d44SBarry Smith      D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1266e69d4d44SBarry Smith      option.
1267e69d4d44SBarry Smith 
1268edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1269edf189efSBarry Smith      is used automatically for a second block.
1270edf189efSBarry Smith 
1271a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
12720971522cSBarry Smith 
1273a541d17aSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners
1274e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
12750971522cSBarry Smith M*/
12760971522cSBarry Smith 
12770971522cSBarry Smith EXTERN_C_BEGIN
12780971522cSBarry Smith #undef __FUNCT__
12790971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
1280*7087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
12810971522cSBarry Smith {
12820971522cSBarry Smith   PetscErrorCode ierr;
12830971522cSBarry Smith   PC_FieldSplit  *jac;
12840971522cSBarry Smith 
12850971522cSBarry Smith   PetscFunctionBegin;
128638f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
12870971522cSBarry Smith   jac->bs        = -1;
12880971522cSBarry Smith   jac->nsplits   = 0;
12893e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1290e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1291c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
129251f519a2SBarry Smith 
12930971522cSBarry Smith   pc->data     = (void*)jac;
12940971522cSBarry Smith 
12950971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1296421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
12970971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
12980971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
12990971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
13000971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
13010971522cSBarry Smith   pc->ops->applyrichardson   = 0;
13020971522cSBarry Smith 
130369a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
130469a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13050971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
13060971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1307704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1308704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
130979416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
131079416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
131151f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
131251f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
13130971522cSBarry Smith   PetscFunctionReturn(0);
13140971522cSBarry Smith }
13150971522cSBarry Smith EXTERN_C_END
13160971522cSBarry Smith 
13170971522cSBarry Smith 
1318a541d17aSBarry Smith /*MC
1319a541d17aSBarry Smith    Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an
1320a541d17aSBarry Smith           overview of these methods.
1321a541d17aSBarry Smith 
1322a541d17aSBarry Smith       Consider the solution to ( A B ) (x_1)  =  (b_1)
1323a541d17aSBarry Smith                                ( C D ) (x_2)     (b_2)
1324a541d17aSBarry Smith 
1325a541d17aSBarry Smith       Important special cases, the Stokes equation: C = B' and D = 0  (A   B) (x_1) = (b_1)
1326a541d17aSBarry Smith                                                                        B'  0) (x_2)   (b_2)
1327a541d17aSBarry Smith 
1328a541d17aSBarry Smith       One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners
1329a541d17aSBarry Smith       for this block system.
1330a541d17aSBarry Smith 
1331a541d17aSBarry Smith       Consider an additional matrix (Ap  Bp)
1332a541d17aSBarry Smith                                     (Cp  Dp) where some or all of the entries may be the same as
1333a541d17aSBarry Smith       in the original matrix (for example Ap == A).
1334a541d17aSBarry Smith 
1335a541d17aSBarry Smith       In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the
1336a541d17aSBarry Smith       approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...)
1337a541d17aSBarry Smith 
1338a541d17aSBarry Smith       Block Jacobi:   x_1 = A^ b_1
1339a541d17aSBarry Smith                       x_2 = D^ b_2
1340a541d17aSBarry Smith 
1341a541d17aSBarry Smith       Lower block Gauss-Seidel:   x_1 = A^ b_1
1342a541d17aSBarry Smith                             x_2 = D^ (b_2 - C x_1)       variant x_2 = D^ (b_2 - Cp x_1)
1343a541d17aSBarry Smith 
1344a541d17aSBarry 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)
1345a541d17aSBarry Smith           Interestingly this form is not actually a symmetric matrix, the symmetric version is
1346a541d17aSBarry Smith                               x_1 = A^(b_1 - B x_2)      variant x_1 = A^(b_1 - Bp x_2)
1347a541d17aSBarry Smith 
13480bc0a719SBarry Smith    Level: intermediate
13490bc0a719SBarry Smith 
1350a541d17aSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT
1351a541d17aSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1352a541d17aSBarry Smith M*/
1353