xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 3c48a1e8da19189ff2402a4e41a2fc082d52c349)
1dba47a55SKris Buschelman 
26356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
3*3c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.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 */
42a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
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     }
145435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 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);
154435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \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;
2086ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = 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);
240435f959eSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
2416ce1633cSBarry Smith       if (stokes) {
2426ce1633cSBarry Smith         IS       zerodiags,rest;
2436ce1633cSBarry Smith         PetscInt nmin,nmax;
2446ce1633cSBarry Smith 
2456ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2466ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2476ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2486ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2496ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
2506ce1633cSBarry Smith         ierr = ISDestroy(zerodiags);CHKERRQ(ierr);
2516ce1633cSBarry Smith         ierr = ISDestroy(rest);CHKERRQ(ierr);
2526ce1633cSBarry Smith       } else {
253d32f9abdSBarry Smith         if (!flg) {
254d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
255d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2566c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2576c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
258d32f9abdSBarry Smith         }
2596c924f48SJed Brown         if (flg || !jac->splitdefined) {
260d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
261db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2626c924f48SJed Brown             char splitname[8];
2636c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
264db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
26579416396SBarry Smith           }
26697bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
267521d7252SBarry Smith         }
26866ffff09SJed Brown       }
2696ce1633cSBarry Smith     }
270edf189efSBarry Smith   } else if (jac->nsplits == 1) {
271edf189efSBarry Smith     if (ilink->is) {
272edf189efSBarry Smith       IS       is2;
273edf189efSBarry Smith       PetscInt nmin,nmax;
274edf189efSBarry Smith 
275edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
276edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
277db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
278edf189efSBarry Smith       ierr = ISDestroy(is2);CHKERRQ(ierr);
279e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
280edf189efSBarry Smith   }
281e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
28269a612a9SBarry Smith   PetscFunctionReturn(0);
28369a612a9SBarry Smith }
28469a612a9SBarry Smith 
28569a612a9SBarry Smith 
28669a612a9SBarry Smith #undef __FUNCT__
28769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
28869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
28969a612a9SBarry Smith {
29069a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
29169a612a9SBarry Smith   PetscErrorCode    ierr;
2925a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
293cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
29469a612a9SBarry Smith   MatStructure      flag = pc->flag;
295ace3abfcSBarry Smith   PetscBool         sorted;
29669a612a9SBarry Smith 
29769a612a9SBarry Smith   PetscFunctionBegin;
29869a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
29997bbdb24SBarry Smith   nsplit = jac->nsplits;
3005a9f2f41SSatish Balay   ilink  = jac->head;
30197bbdb24SBarry Smith 
30297bbdb24SBarry Smith   /* get the matrices for each split */
303704ba839SBarry Smith   if (!jac->issetup) {
3041b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
30597bbdb24SBarry Smith 
306704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
307704ba839SBarry Smith 
308704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
30951f519a2SBarry Smith     bs     = jac->bs;
31097bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
311cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3121b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3131b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3141b9fc7fcSBarry Smith       if (jac->defaultsplit) {
315704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
316704ba839SBarry Smith       } else if (!ilink->is) {
317ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3185a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3195a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3201b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3211b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3221b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
32397bbdb24SBarry Smith             }
32497bbdb24SBarry Smith           }
325d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
326ccb205f8SBarry Smith         } else {
327704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
328ccb205f8SBarry Smith         }
3293e197d65SBarry Smith       }
330704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
331e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
332704ba839SBarry Smith       ilink = ilink->next;
3331b9fc7fcSBarry Smith     }
3341b9fc7fcSBarry Smith   }
3351b9fc7fcSBarry Smith 
336704ba839SBarry Smith   ilink  = jac->head;
33797bbdb24SBarry Smith   if (!jac->pmat) {
338cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
339cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3404aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
341704ba839SBarry Smith       ilink = ilink->next;
342cf502942SBarry Smith     }
34397bbdb24SBarry Smith   } else {
344cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3454aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
346704ba839SBarry Smith       ilink = ilink->next;
347cf502942SBarry Smith     }
34897bbdb24SBarry Smith   }
349519d70e2SJed Brown   if (jac->realdiagonal) {
350519d70e2SJed Brown     ilink = jac->head;
351519d70e2SJed Brown     if (!jac->mat) {
352519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
353519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
354519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
355519d70e2SJed Brown         ilink = ilink->next;
356519d70e2SJed Brown       }
357519d70e2SJed Brown     } else {
358519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
359519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
360519d70e2SJed Brown         ilink = ilink->next;
361519d70e2SJed Brown       }
362519d70e2SJed Brown     }
363519d70e2SJed Brown   } else {
364519d70e2SJed Brown     jac->mat = jac->pmat;
365519d70e2SJed Brown   }
36697bbdb24SBarry Smith 
3676c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
36868dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
36968dd23aaSBarry Smith     ilink  = jac->head;
37068dd23aaSBarry Smith     if (!jac->Afield) {
37168dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
37268dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3734aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
37468dd23aaSBarry Smith         ilink = ilink->next;
37568dd23aaSBarry Smith       }
37668dd23aaSBarry Smith     } else {
37768dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3784aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
37968dd23aaSBarry Smith         ilink = ilink->next;
38068dd23aaSBarry Smith       }
38168dd23aaSBarry Smith     }
38268dd23aaSBarry Smith   }
38368dd23aaSBarry Smith 
3843b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
3853b224e63SBarry Smith     IS       ccis;
3864aa3045dSJed Brown     PetscInt rstart,rend;
387e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
38868dd23aaSBarry Smith 
389e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
390e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
391e6cab6aaSJed Brown 
3923b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
3933b224e63SBarry Smith     if (jac->schur) {
3943b224e63SBarry Smith       ilink = jac->head;
3954aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3964aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
3973b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3983b224e63SBarry Smith       ilink = ilink->next;
3994aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4004aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
4013b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
402519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
403084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4043b224e63SBarry Smith 
4053b224e63SBarry Smith      } else {
4061cee3971SBarry Smith       KSP ksp;
4076c924f48SJed Brown       char schurprefix[256];
4083b224e63SBarry Smith 
409a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4103b224e63SBarry Smith       ilink = jac->head;
4114aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4124aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
4133b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
4143b224e63SBarry Smith       ilink = ilink->next;
4154aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4164aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
4173b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
418084e4875SJed Brown       /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */
419519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
420a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4211cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
422e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
423a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
424a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
4253b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4263b224e63SBarry Smith 
4273b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4289005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4291cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
430084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
431084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
432084e4875SJed Brown         PC pc;
433cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
434084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
435084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
436e69d4d44SBarry Smith       }
4376c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4386c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4393b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4403b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4413b224e63SBarry Smith 
4423b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4433b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4443b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4453b224e63SBarry Smith       ilink = jac->head;
4463b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4473b224e63SBarry Smith       ilink = ilink->next;
4483b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4493b224e63SBarry Smith     }
4503b224e63SBarry Smith   } else {
45197bbdb24SBarry Smith     /* set up the individual PCs */
45297bbdb24SBarry Smith     i    = 0;
4535a9f2f41SSatish Balay     ilink = jac->head;
4545a9f2f41SSatish Balay     while (ilink) {
455519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
4563b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4575a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
4585a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
45997bbdb24SBarry Smith       i++;
4605a9f2f41SSatish Balay       ilink = ilink->next;
4610971522cSBarry Smith     }
46297bbdb24SBarry Smith 
46397bbdb24SBarry Smith     /* create work vectors for each split */
4641b9fc7fcSBarry Smith     if (!jac->x) {
46597bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
4665a9f2f41SSatish Balay       ilink = jac->head;
46797bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
468906ed7ccSBarry Smith         Vec *vl,*vr;
469906ed7ccSBarry Smith 
470906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
471906ed7ccSBarry Smith         ilink->x  = *vr;
472906ed7ccSBarry Smith         ilink->y  = *vl;
473906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
474906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
4755a9f2f41SSatish Balay         jac->x[i] = ilink->x;
4765a9f2f41SSatish Balay         jac->y[i] = ilink->y;
4775a9f2f41SSatish Balay         ilink     = ilink->next;
47897bbdb24SBarry Smith       }
4793b224e63SBarry Smith     }
4803b224e63SBarry Smith   }
4813b224e63SBarry Smith 
4823b224e63SBarry Smith 
483d0f46423SBarry Smith   if (!jac->head->sctx) {
4843b224e63SBarry Smith     Vec xtmp;
4853b224e63SBarry Smith 
48679416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
4871b9fc7fcSBarry Smith 
4885a9f2f41SSatish Balay     ilink = jac->head;
4891b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
4901b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
491704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
4925a9f2f41SSatish Balay       ilink = ilink->next;
4931b9fc7fcSBarry Smith     }
4941b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
4951b9fc7fcSBarry Smith   }
4960971522cSBarry Smith   PetscFunctionReturn(0);
4970971522cSBarry Smith }
4980971522cSBarry Smith 
4995a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
500ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
501ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5025a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
503ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
504ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
50579416396SBarry Smith 
5060971522cSBarry Smith #undef __FUNCT__
5073b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5083b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5093b224e63SBarry Smith {
5103b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5113b224e63SBarry Smith   PetscErrorCode    ierr;
5123b224e63SBarry Smith   KSP               ksp;
5133b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5143b224e63SBarry Smith 
5153b224e63SBarry Smith   PetscFunctionBegin;
5163b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5173b224e63SBarry Smith 
518c5d2311dSJed Brown   switch (jac->schurfactorization) {
519c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
520a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
521c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
522c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
523c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
524c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
525c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
526c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
527c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
528c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
529c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
531c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
532c5d2311dSJed Brown       break;
533c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
534a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
535c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
536c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
537c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
538c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
539c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
540c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
541c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
542c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
543c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
544c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
545c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
546c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
547c5d2311dSJed Brown       break;
548c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
549a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
550c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
551c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
552c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
553c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
554c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
555c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
556c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
557c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
558c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
559c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
560c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
561c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
562c5d2311dSJed Brown       break;
563c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
564a04f6461SBarry Smith       /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1], an exact solve if applied exactly, needs one extra solve with A */
5653b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5663b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5673b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5683b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
5693b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
5703b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5713b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5723b224e63SBarry Smith 
5733b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
5743b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5753b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5763b224e63SBarry Smith 
5773b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
5783b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
5793b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5803b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5813b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
582c5d2311dSJed Brown   }
5833b224e63SBarry Smith   PetscFunctionReturn(0);
5843b224e63SBarry Smith }
5853b224e63SBarry Smith 
5863b224e63SBarry Smith #undef __FUNCT__
5870971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
5880971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
5890971522cSBarry Smith {
5900971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5910971522cSBarry Smith   PetscErrorCode    ierr;
5925a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
593af93645bSJed Brown   PetscInt          cnt;
5940971522cSBarry Smith 
5950971522cSBarry Smith   PetscFunctionBegin;
59651f519a2SBarry Smith   CHKMEMQ;
59751f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
59851f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
59951f519a2SBarry Smith 
60079416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6011b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6020971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6035a9f2f41SSatish Balay       while (ilink) {
6045a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6055a9f2f41SSatish Balay         ilink = ilink->next;
6060971522cSBarry Smith       }
6070971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6081b9fc7fcSBarry Smith     } else {
609efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6105a9f2f41SSatish Balay       while (ilink) {
6115a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6125a9f2f41SSatish Balay         ilink = ilink->next;
6131b9fc7fcSBarry Smith       }
6141b9fc7fcSBarry Smith     }
61516913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
61679416396SBarry Smith     if (!jac->w1) {
61779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
61879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
61979416396SBarry Smith     }
620efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6215a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6223e197d65SBarry Smith     cnt = 1;
6235a9f2f41SSatish Balay     while (ilink->next) {
6245a9f2f41SSatish Balay       ilink = ilink->next;
6253e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6263e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6273e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6283e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6293e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6303e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6313e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6323e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6333e197d65SBarry Smith     }
63451f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
63511755939SBarry Smith       cnt -= 2;
63651f519a2SBarry Smith       while (ilink->previous) {
63751f519a2SBarry Smith         ilink = ilink->previous;
63811755939SBarry Smith         /* compute the residual only over the part of the vector needed */
63911755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
64011755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
64111755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64211755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64311755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
64411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64651f519a2SBarry Smith       }
64711755939SBarry Smith     }
64865e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
64951f519a2SBarry Smith   CHKMEMQ;
6500971522cSBarry Smith   PetscFunctionReturn(0);
6510971522cSBarry Smith }
6520971522cSBarry Smith 
653421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
654ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
655ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
656421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
657ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
658ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
659421e10b8SBarry Smith 
660421e10b8SBarry Smith #undef __FUNCT__
6618c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
662421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
663421e10b8SBarry Smith {
664421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
665421e10b8SBarry Smith   PetscErrorCode    ierr;
666421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
667421e10b8SBarry Smith 
668421e10b8SBarry Smith   PetscFunctionBegin;
669421e10b8SBarry Smith   CHKMEMQ;
670421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
671421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
672421e10b8SBarry Smith 
673421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
674421e10b8SBarry Smith     if (jac->defaultsplit) {
675421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
676421e10b8SBarry Smith       while (ilink) {
677421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
678421e10b8SBarry Smith 	ilink = ilink->next;
679421e10b8SBarry Smith       }
680421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
681421e10b8SBarry Smith     } else {
682421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
683421e10b8SBarry Smith       while (ilink) {
684421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
685421e10b8SBarry Smith 	ilink = ilink->next;
686421e10b8SBarry Smith       }
687421e10b8SBarry Smith     }
688421e10b8SBarry Smith   } else {
689421e10b8SBarry Smith     if (!jac->w1) {
690421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
691421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
692421e10b8SBarry Smith     }
693421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
694421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
695421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
696421e10b8SBarry Smith       while (ilink->next) {
697421e10b8SBarry Smith         ilink = ilink->next;
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       while (ilink->previous) {
703421e10b8SBarry Smith         ilink = ilink->previous;
7049989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
705421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
706421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
707421e10b8SBarry Smith       }
708421e10b8SBarry Smith     } else {
709421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
710421e10b8SBarry Smith 	ilink = ilink->next;
711421e10b8SBarry Smith       }
712421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
713421e10b8SBarry Smith       while (ilink->previous) {
714421e10b8SBarry Smith 	ilink = ilink->previous;
7159989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
716421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
717421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
718421e10b8SBarry Smith       }
719421e10b8SBarry Smith     }
720421e10b8SBarry Smith   }
721421e10b8SBarry Smith   CHKMEMQ;
722421e10b8SBarry Smith   PetscFunctionReturn(0);
723421e10b8SBarry Smith }
724421e10b8SBarry Smith 
7250971522cSBarry Smith #undef __FUNCT__
726574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
727574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
7280971522cSBarry Smith {
7290971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7300971522cSBarry Smith   PetscErrorCode    ierr;
7315a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7320971522cSBarry Smith 
7330971522cSBarry Smith   PetscFunctionBegin;
7345a9f2f41SSatish Balay   while (ilink) {
735574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
7365a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
7375a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
7385a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
739704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
7405a9f2f41SSatish Balay     next = ilink->next;
7415a9f2f41SSatish Balay     ilink = next;
7420971522cSBarry Smith   }
74305b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
744574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
745574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
746574deadeSBarry Smith   } else if (jac->mat) {
747574deadeSBarry Smith     jac->mat = PETSC_NULL;
748574deadeSBarry Smith   }
74997bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
75068dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
75179416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
75279416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
7533b224e63SBarry Smith   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
754084e4875SJed Brown   if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
755d0f46423SBarry Smith   if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);}
7563b224e63SBarry Smith   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
7573b224e63SBarry Smith   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
758574deadeSBarry Smith   PetscFunctionReturn(0);
759574deadeSBarry Smith }
760574deadeSBarry Smith 
761574deadeSBarry Smith #undef __FUNCT__
762574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
763574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
764574deadeSBarry Smith {
765574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
766574deadeSBarry Smith   PetscErrorCode    ierr;
767574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
768574deadeSBarry Smith 
769574deadeSBarry Smith   PetscFunctionBegin;
770574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
771574deadeSBarry Smith   while (ilink) {
772574deadeSBarry Smith     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
773574deadeSBarry Smith     next = ilink->next;
774574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
775574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
776574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
777574deadeSBarry Smith     ilink = next;
778574deadeSBarry Smith   }
779574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
780c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
7810971522cSBarry Smith   PetscFunctionReturn(0);
7820971522cSBarry Smith }
7830971522cSBarry Smith 
7840971522cSBarry Smith #undef __FUNCT__
7850971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
7860971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
7870971522cSBarry Smith {
7881b9fc7fcSBarry Smith   PetscErrorCode  ierr;
7896c924f48SJed Brown   PetscInt        bs;
790bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
7919dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
7923b224e63SBarry Smith   PCCompositeType ctype;
7931b9fc7fcSBarry Smith 
7940971522cSBarry Smith   PetscFunctionBegin;
7951b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
796acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
79751f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
79851f519a2SBarry Smith   if (flg) {
79951f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
80051f519a2SBarry Smith   }
801704ba839SBarry Smith 
802435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
803c0adfefeSBarry Smith   if (stokes) {
804c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
805c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
806c0adfefeSBarry Smith   }
807c0adfefeSBarry Smith 
8083b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8093b224e63SBarry Smith   if (flg) {
8103b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8113b224e63SBarry Smith   }
812d32f9abdSBarry Smith 
813c30613efSMatthew Knepley   /* Only setup fields once */
814c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
815d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
816d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8176c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8186c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
819d32f9abdSBarry Smith   }
820c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
821c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
822084e4875SJed 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);
823c5d2311dSJed Brown   }
8241b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8250971522cSBarry Smith   PetscFunctionReturn(0);
8260971522cSBarry Smith }
8270971522cSBarry Smith 
8280971522cSBarry Smith /*------------------------------------------------------------------------------------*/
8290971522cSBarry Smith 
8300971522cSBarry Smith EXTERN_C_BEGIN
8310971522cSBarry Smith #undef __FUNCT__
8320971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
8337087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
8340971522cSBarry Smith {
83597bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8360971522cSBarry Smith   PetscErrorCode    ierr;
8375a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
83869a612a9SBarry Smith   char              prefix[128];
83951f519a2SBarry Smith   PetscInt          i;
8400971522cSBarry Smith 
8410971522cSBarry Smith   PetscFunctionBegin;
8426c924f48SJed Brown   if (jac->splitdefined) {
8436c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8446c924f48SJed Brown     PetscFunctionReturn(0);
8456c924f48SJed Brown   }
84651f519a2SBarry Smith   for (i=0; i<n; i++) {
847e32f2f54SBarry 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);
848e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
84951f519a2SBarry Smith   }
850704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
851a04f6461SBarry Smith   if (splitname) {
852db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
853a04f6461SBarry Smith   } else {
854a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
855a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
856a04f6461SBarry Smith   }
857704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
8585a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
8595a9f2f41SSatish Balay   ilink->nfields = n;
8605a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
8617adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
8621cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
8635a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
8649005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
86569a612a9SBarry Smith 
866a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
8675a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
8680971522cSBarry Smith 
8690971522cSBarry Smith   if (!next) {
8705a9f2f41SSatish Balay     jac->head       = ilink;
87151f519a2SBarry Smith     ilink->previous = PETSC_NULL;
8720971522cSBarry Smith   } else {
8730971522cSBarry Smith     while (next->next) {
8740971522cSBarry Smith       next = next->next;
8750971522cSBarry Smith     }
8765a9f2f41SSatish Balay     next->next      = ilink;
87751f519a2SBarry Smith     ilink->previous = next;
8780971522cSBarry Smith   }
8790971522cSBarry Smith   jac->nsplits++;
8800971522cSBarry Smith   PetscFunctionReturn(0);
8810971522cSBarry Smith }
8820971522cSBarry Smith EXTERN_C_END
8830971522cSBarry Smith 
884e69d4d44SBarry Smith EXTERN_C_BEGIN
885e69d4d44SBarry Smith #undef __FUNCT__
886e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
8877087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
888e69d4d44SBarry Smith {
889e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
890e69d4d44SBarry Smith   PetscErrorCode ierr;
891e69d4d44SBarry Smith 
892e69d4d44SBarry Smith   PetscFunctionBegin;
893519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
894e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
895e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
89613e0d083SBarry Smith   if (n) *n = jac->nsplits;
897e69d4d44SBarry Smith   PetscFunctionReturn(0);
898e69d4d44SBarry Smith }
899e69d4d44SBarry Smith EXTERN_C_END
9000971522cSBarry Smith 
9010971522cSBarry Smith EXTERN_C_BEGIN
9020971522cSBarry Smith #undef __FUNCT__
90369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9047087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9050971522cSBarry Smith {
9060971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9070971522cSBarry Smith   PetscErrorCode    ierr;
9080971522cSBarry Smith   PetscInt          cnt = 0;
9095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9100971522cSBarry Smith 
9110971522cSBarry Smith   PetscFunctionBegin;
91269a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
9135a9f2f41SSatish Balay   while (ilink) {
9145a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9155a9f2f41SSatish Balay     ilink = ilink->next;
9160971522cSBarry Smith   }
917e32f2f54SBarry 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);
91813e0d083SBarry Smith   if (n) *n = jac->nsplits;
9190971522cSBarry Smith   PetscFunctionReturn(0);
9200971522cSBarry Smith }
9210971522cSBarry Smith EXTERN_C_END
9220971522cSBarry Smith 
923704ba839SBarry Smith EXTERN_C_BEGIN
924704ba839SBarry Smith #undef __FUNCT__
925704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
9267087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
927704ba839SBarry Smith {
928704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
929704ba839SBarry Smith   PetscErrorCode    ierr;
930704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
931704ba839SBarry Smith   char              prefix[128];
932704ba839SBarry Smith 
933704ba839SBarry Smith   PetscFunctionBegin;
9346c924f48SJed Brown   if (jac->splitdefined) {
9356c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9366c924f48SJed Brown     PetscFunctionReturn(0);
9376c924f48SJed Brown   }
93816913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
939a04f6461SBarry Smith   if (splitname) {
940db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
941a04f6461SBarry Smith   } else {
942a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
943a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
944a04f6461SBarry Smith   }
9451ab39975SBarry Smith   ilink->is      = is;
9461ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
947704ba839SBarry Smith   ilink->next    = PETSC_NULL;
948704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9491cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
950704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9519005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
952704ba839SBarry Smith 
953a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
954704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
955704ba839SBarry Smith 
956704ba839SBarry Smith   if (!next) {
957704ba839SBarry Smith     jac->head       = ilink;
958704ba839SBarry Smith     ilink->previous = PETSC_NULL;
959704ba839SBarry Smith   } else {
960704ba839SBarry Smith     while (next->next) {
961704ba839SBarry Smith       next = next->next;
962704ba839SBarry Smith     }
963704ba839SBarry Smith     next->next      = ilink;
964704ba839SBarry Smith     ilink->previous = next;
965704ba839SBarry Smith   }
966704ba839SBarry Smith   jac->nsplits++;
967704ba839SBarry Smith 
968704ba839SBarry Smith   PetscFunctionReturn(0);
969704ba839SBarry Smith }
970704ba839SBarry Smith EXTERN_C_END
971704ba839SBarry Smith 
9720971522cSBarry Smith #undef __FUNCT__
9730971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
9740971522cSBarry Smith /*@
9750971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
9760971522cSBarry Smith 
977ad4df100SBarry Smith     Logically Collective on PC
9780971522cSBarry Smith 
9790971522cSBarry Smith     Input Parameters:
9800971522cSBarry Smith +   pc  - the preconditioner context
981a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
9820971522cSBarry Smith .   n - the number of fields in this split
983db4c96c1SJed Brown -   fields - the fields in this split
9840971522cSBarry Smith 
9850971522cSBarry Smith     Level: intermediate
9860971522cSBarry Smith 
987d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
988d32f9abdSBarry Smith 
989d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
990d32f9abdSBarry 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
991d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
992d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
993d32f9abdSBarry Smith 
994db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
995db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
996db4c96c1SJed Brown 
997d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
9980971522cSBarry Smith 
9990971522cSBarry Smith @*/
10007087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10010971522cSBarry Smith {
10024ac538c5SBarry Smith   PetscErrorCode ierr;
10030971522cSBarry Smith 
10040971522cSBarry Smith   PetscFunctionBegin;
10050700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1006db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1007db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1008db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10094ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10100971522cSBarry Smith   PetscFunctionReturn(0);
10110971522cSBarry Smith }
10120971522cSBarry Smith 
10130971522cSBarry Smith #undef __FUNCT__
1014704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1015704ba839SBarry Smith /*@
1016704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1017704ba839SBarry Smith 
1018ad4df100SBarry Smith     Logically Collective on PC
1019704ba839SBarry Smith 
1020704ba839SBarry Smith     Input Parameters:
1021704ba839SBarry Smith +   pc  - the preconditioner context
1022a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1023db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1024704ba839SBarry Smith 
1025d32f9abdSBarry Smith 
1026a6ffb8dbSJed Brown     Notes:
1027a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1028a6ffb8dbSJed Brown 
1029db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1030db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1031d32f9abdSBarry Smith 
1032704ba839SBarry Smith     Level: intermediate
1033704ba839SBarry Smith 
1034704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1035704ba839SBarry Smith 
1036704ba839SBarry Smith @*/
10377087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1038704ba839SBarry Smith {
10394ac538c5SBarry Smith   PetscErrorCode ierr;
1040704ba839SBarry Smith 
1041704ba839SBarry Smith   PetscFunctionBegin;
10420700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1043db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1044db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
10454ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1046704ba839SBarry Smith   PetscFunctionReturn(0);
1047704ba839SBarry Smith }
1048704ba839SBarry Smith 
1049704ba839SBarry Smith #undef __FUNCT__
105051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
105151f519a2SBarry Smith /*@
105251f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
105351f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
105451f519a2SBarry Smith 
1055ad4df100SBarry Smith     Logically Collective on PC
105651f519a2SBarry Smith 
105751f519a2SBarry Smith     Input Parameters:
105851f519a2SBarry Smith +   pc  - the preconditioner context
105951f519a2SBarry Smith -   bs - the block size
106051f519a2SBarry Smith 
106151f519a2SBarry Smith     Level: intermediate
106251f519a2SBarry Smith 
106351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
106451f519a2SBarry Smith 
106551f519a2SBarry Smith @*/
10667087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
106751f519a2SBarry Smith {
10684ac538c5SBarry Smith   PetscErrorCode ierr;
106951f519a2SBarry Smith 
107051f519a2SBarry Smith   PetscFunctionBegin;
10710700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1072c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
10734ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
107451f519a2SBarry Smith   PetscFunctionReturn(0);
107551f519a2SBarry Smith }
107651f519a2SBarry Smith 
107751f519a2SBarry Smith #undef __FUNCT__
107869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
10790971522cSBarry Smith /*@C
108069a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
10810971522cSBarry Smith 
108269a612a9SBarry Smith    Collective on KSP
10830971522cSBarry Smith 
10840971522cSBarry Smith    Input Parameter:
10850971522cSBarry Smith .  pc - the preconditioner context
10860971522cSBarry Smith 
10870971522cSBarry Smith    Output Parameters:
108813e0d083SBarry Smith +  n - the number of splits
108969a612a9SBarry Smith -  pc - the array of KSP contexts
10900971522cSBarry Smith 
10910971522cSBarry Smith    Note:
1092d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1093d32f9abdSBarry Smith    (not the KSP just the array that contains them).
10940971522cSBarry Smith 
109569a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
10960971522cSBarry Smith 
10970971522cSBarry Smith    Level: advanced
10980971522cSBarry Smith 
10990971522cSBarry Smith .seealso: PCFIELDSPLIT
11000971522cSBarry Smith @*/
11017087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
11020971522cSBarry Smith {
11034ac538c5SBarry Smith   PetscErrorCode ierr;
11040971522cSBarry Smith 
11050971522cSBarry Smith   PetscFunctionBegin;
11060700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
110713e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
11084ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
11090971522cSBarry Smith   PetscFunctionReturn(0);
11100971522cSBarry Smith }
11110971522cSBarry Smith 
1112e69d4d44SBarry Smith #undef __FUNCT__
1113e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1114e69d4d44SBarry Smith /*@
1115e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1116a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1117e69d4d44SBarry Smith 
1118e69d4d44SBarry Smith     Collective on PC
1119e69d4d44SBarry Smith 
1120e69d4d44SBarry Smith     Input Parameters:
1121e69d4d44SBarry Smith +   pc  - the preconditioner context
1122084e4875SJed Brown .   ptype - which matrix to use for preconditioning the Schur complement
1123084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1124084e4875SJed Brown 
1125084e4875SJed Brown     Notes:
1126a04f6461SBarry Smith     The default is to use the block on the diagonal of the preconditioning matrix.  This is A11, in the (1,1) position.
1127084e4875SJed Brown     There are currently no preconditioners that work directly with the Schur complement so setting
1128084e4875SJed Brown     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
1129e69d4d44SBarry Smith 
1130e69d4d44SBarry Smith     Options Database:
1131084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1132e69d4d44SBarry Smith 
1133e69d4d44SBarry Smith     Level: intermediate
1134e69d4d44SBarry Smith 
1135084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1136e69d4d44SBarry Smith 
1137e69d4d44SBarry Smith @*/
11387087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1139e69d4d44SBarry Smith {
11404ac538c5SBarry Smith   PetscErrorCode ierr;
1141e69d4d44SBarry Smith 
1142e69d4d44SBarry Smith   PetscFunctionBegin;
11430700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11444ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1145e69d4d44SBarry Smith   PetscFunctionReturn(0);
1146e69d4d44SBarry Smith }
1147e69d4d44SBarry Smith 
1148e69d4d44SBarry Smith EXTERN_C_BEGIN
1149e69d4d44SBarry Smith #undef __FUNCT__
1150e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
11517087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1152e69d4d44SBarry Smith {
1153e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1154084e4875SJed Brown   PetscErrorCode  ierr;
1155e69d4d44SBarry Smith 
1156e69d4d44SBarry Smith   PetscFunctionBegin;
1157084e4875SJed Brown   jac->schurpre = ptype;
1158084e4875SJed Brown   if (pre) {
1159084e4875SJed Brown     if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
1160084e4875SJed Brown     jac->schur_user = pre;
1161084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1162084e4875SJed Brown   }
1163e69d4d44SBarry Smith   PetscFunctionReturn(0);
1164e69d4d44SBarry Smith }
1165e69d4d44SBarry Smith EXTERN_C_END
1166e69d4d44SBarry Smith 
116730ad9308SMatthew Knepley #undef __FUNCT__
116830ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
116930ad9308SMatthew Knepley /*@C
11708c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
117130ad9308SMatthew Knepley 
117230ad9308SMatthew Knepley    Collective on KSP
117330ad9308SMatthew Knepley 
117430ad9308SMatthew Knepley    Input Parameter:
117530ad9308SMatthew Knepley .  pc - the preconditioner context
117630ad9308SMatthew Knepley 
117730ad9308SMatthew Knepley    Output Parameters:
1178a04f6461SBarry Smith +  A00 - the (0,0) block
1179a04f6461SBarry Smith .  A01 - the (0,1) block
1180a04f6461SBarry Smith .  A10 - the (1,0) block
1181a04f6461SBarry Smith -  A11 - the (1,1) block
118230ad9308SMatthew Knepley 
118330ad9308SMatthew Knepley    Level: advanced
118430ad9308SMatthew Knepley 
118530ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
118630ad9308SMatthew Knepley @*/
1187a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
118830ad9308SMatthew Knepley {
118930ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
119030ad9308SMatthew Knepley 
119130ad9308SMatthew Knepley   PetscFunctionBegin;
11920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1193c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1194a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1195a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1196a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1197a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
119830ad9308SMatthew Knepley   PetscFunctionReturn(0);
119930ad9308SMatthew Knepley }
120030ad9308SMatthew Knepley 
120179416396SBarry Smith EXTERN_C_BEGIN
120279416396SBarry Smith #undef __FUNCT__
120379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
12047087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
120579416396SBarry Smith {
120679416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1207e69d4d44SBarry Smith   PetscErrorCode ierr;
120879416396SBarry Smith 
120979416396SBarry Smith   PetscFunctionBegin;
121079416396SBarry Smith   jac->type = type;
12113b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
12123b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
12133b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1214e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1215e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1216e69d4d44SBarry Smith 
12173b224e63SBarry Smith   } else {
12183b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
12193b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1220e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
12219e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
12223b224e63SBarry Smith   }
122379416396SBarry Smith   PetscFunctionReturn(0);
122479416396SBarry Smith }
122579416396SBarry Smith EXTERN_C_END
122679416396SBarry Smith 
122751f519a2SBarry Smith EXTERN_C_BEGIN
122851f519a2SBarry Smith #undef __FUNCT__
122951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
12307087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
123151f519a2SBarry Smith {
123251f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
123351f519a2SBarry Smith 
123451f519a2SBarry Smith   PetscFunctionBegin;
123565e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
123665e19b50SBarry 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);
123751f519a2SBarry Smith   jac->bs = bs;
123851f519a2SBarry Smith   PetscFunctionReturn(0);
123951f519a2SBarry Smith }
124051f519a2SBarry Smith EXTERN_C_END
124151f519a2SBarry Smith 
124279416396SBarry Smith #undef __FUNCT__
124379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1244bc08b0f1SBarry Smith /*@
124579416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
124679416396SBarry Smith 
124779416396SBarry Smith    Collective on PC
124879416396SBarry Smith 
124979416396SBarry Smith    Input Parameter:
125079416396SBarry Smith .  pc - the preconditioner context
125181540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
125279416396SBarry Smith 
125379416396SBarry Smith    Options Database Key:
1254a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
125579416396SBarry Smith 
125679416396SBarry Smith    Level: Developer
125779416396SBarry Smith 
125879416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
125979416396SBarry Smith 
126079416396SBarry Smith .seealso: PCCompositeSetType()
126179416396SBarry Smith 
126279416396SBarry Smith @*/
12637087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
126479416396SBarry Smith {
12654ac538c5SBarry Smith   PetscErrorCode ierr;
126679416396SBarry Smith 
126779416396SBarry Smith   PetscFunctionBegin;
12680700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12694ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
127079416396SBarry Smith   PetscFunctionReturn(0);
127179416396SBarry Smith }
127279416396SBarry Smith 
12730971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
12740971522cSBarry Smith /*MC
1275a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1276a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
12770971522cSBarry Smith 
1278edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1279edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
12800971522cSBarry Smith 
1281a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
128269a612a9SBarry Smith          and set the options directly on the resulting KSP object
12830971522cSBarry Smith 
12840971522cSBarry Smith    Level: intermediate
12850971522cSBarry Smith 
128679416396SBarry Smith    Options Database Keys:
128781540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
128881540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
128981540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
129081540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
129181540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1292e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1293435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
129479416396SBarry Smith 
1295edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
12963b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
12973b224e63SBarry Smith 
1298d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1299d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1300d32f9abdSBarry Smith 
1301d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1302d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1303d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1304d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1305d32f9abdSBarry Smith 
1306a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1307a04f6461SBarry Smith                                                     ( A10 A11 )
1308e69d4d44SBarry Smith      the preconditioner is
1309a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1310a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1311a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1312a04f6461SBarry Smith      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given in providing the SECOND split or 1 if not give).
1313a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1314edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1315a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
13167e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1317e69d4d44SBarry Smith 
1318edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1319edf189efSBarry Smith      is used automatically for a second block.
1320edf189efSBarry Smith 
1321a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
13220971522cSBarry Smith 
13237e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1324e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
13250971522cSBarry Smith M*/
13260971522cSBarry Smith 
13270971522cSBarry Smith EXTERN_C_BEGIN
13280971522cSBarry Smith #undef __FUNCT__
13290971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
13307087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
13310971522cSBarry Smith {
13320971522cSBarry Smith   PetscErrorCode ierr;
13330971522cSBarry Smith   PC_FieldSplit  *jac;
13340971522cSBarry Smith 
13350971522cSBarry Smith   PetscFunctionBegin;
133638f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
13370971522cSBarry Smith   jac->bs        = -1;
13380971522cSBarry Smith   jac->nsplits   = 0;
13393e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1340e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1341c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
134251f519a2SBarry Smith 
13430971522cSBarry Smith   pc->data     = (void*)jac;
13440971522cSBarry Smith 
13450971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1346421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
13470971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1348574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
13490971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
13500971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
13510971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
13520971522cSBarry Smith   pc->ops->applyrichardson   = 0;
13530971522cSBarry Smith 
135469a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
135569a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13560971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
13570971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1358704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1359704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
136079416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
136179416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
136251f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
136351f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
13640971522cSBarry Smith   PetscFunctionReturn(0);
13650971522cSBarry Smith }
13660971522cSBarry Smith EXTERN_C_END
13670971522cSBarry Smith 
13680971522cSBarry Smith 
1369a541d17aSBarry Smith 
1370