xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 5d4804778dfa45ea00bcad8aabdb3eca36eaeac3)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed 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);
226fcfd50ebSBarry Smith           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);
250fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
251fcfd50ebSBarry 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);
278fcfd50ebSBarry 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 #undef __FUNCT__
28669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
28769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
28869a612a9SBarry Smith {
28969a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
29069a612a9SBarry Smith   PetscErrorCode    ierr;
2915a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
292cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
29369a612a9SBarry Smith   MatStructure      flag = pc->flag;
294ace3abfcSBarry Smith   PetscBool         sorted;
29569a612a9SBarry Smith 
29669a612a9SBarry Smith   PetscFunctionBegin;
29769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
29897bbdb24SBarry Smith   nsplit = jac->nsplits;
2995a9f2f41SSatish Balay   ilink  = jac->head;
30097bbdb24SBarry Smith 
30197bbdb24SBarry Smith   /* get the matrices for each split */
302704ba839SBarry Smith   if (!jac->issetup) {
3031b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
30497bbdb24SBarry Smith 
305704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
306704ba839SBarry Smith 
307704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
30851f519a2SBarry Smith     bs     = jac->bs;
30997bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
310cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3111b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3121b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3131b9fc7fcSBarry Smith       if (jac->defaultsplit) {
314704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
315704ba839SBarry Smith       } else if (!ilink->is) {
316ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3175a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3185a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3191b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3201b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3211b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
32297bbdb24SBarry Smith             }
32397bbdb24SBarry Smith           }
324d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
325ccb205f8SBarry Smith         } else {
326704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
327ccb205f8SBarry Smith         }
3283e197d65SBarry Smith       }
329704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
330e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
331704ba839SBarry Smith       ilink = ilink->next;
3321b9fc7fcSBarry Smith     }
3331b9fc7fcSBarry Smith   }
3341b9fc7fcSBarry Smith 
335704ba839SBarry Smith   ilink  = jac->head;
33697bbdb24SBarry Smith   if (!jac->pmat) {
337cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
338cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3394aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
340704ba839SBarry Smith       ilink = ilink->next;
341cf502942SBarry Smith     }
34297bbdb24SBarry Smith   } else {
343cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3444aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
345704ba839SBarry Smith       ilink = ilink->next;
346cf502942SBarry Smith     }
34797bbdb24SBarry Smith   }
348519d70e2SJed Brown   if (jac->realdiagonal) {
349519d70e2SJed Brown     ilink = jac->head;
350519d70e2SJed Brown     if (!jac->mat) {
351519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
352519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
353519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
354519d70e2SJed Brown         ilink = ilink->next;
355519d70e2SJed Brown       }
356519d70e2SJed Brown     } else {
357519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
358966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
359519d70e2SJed Brown         ilink = ilink->next;
360519d70e2SJed Brown       }
361519d70e2SJed Brown     }
362519d70e2SJed Brown   } else {
363519d70e2SJed Brown     jac->mat = jac->pmat;
364519d70e2SJed Brown   }
36597bbdb24SBarry Smith 
3666c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
36768dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
36868dd23aaSBarry Smith     ilink  = jac->head;
36968dd23aaSBarry Smith     if (!jac->Afield) {
37068dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
37168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3724aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
37368dd23aaSBarry Smith         ilink = ilink->next;
37468dd23aaSBarry Smith       }
37568dd23aaSBarry Smith     } else {
37668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3774aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
37868dd23aaSBarry Smith         ilink = ilink->next;
37968dd23aaSBarry Smith       }
38068dd23aaSBarry Smith     }
38168dd23aaSBarry Smith   }
38268dd23aaSBarry Smith 
3833b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
3843b224e63SBarry Smith     IS       ccis;
3854aa3045dSJed Brown     PetscInt rstart,rend;
386e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
38768dd23aaSBarry Smith 
388e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
389e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
390e6cab6aaSJed Brown 
3913b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
3923b224e63SBarry Smith     if (jac->schur) {
3933b224e63SBarry Smith       ilink = jac->head;
3944aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3954aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
396fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
3973b224e63SBarry Smith       ilink = ilink->next;
3984aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3994aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
400fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
401519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
402084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4033b224e63SBarry Smith 
4043b224e63SBarry Smith      } else {
4051cee3971SBarry Smith       KSP ksp;
4066c924f48SJed Brown       char schurprefix[256];
4073b224e63SBarry Smith 
408a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4093b224e63SBarry Smith       ilink = jac->head;
4104aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4114aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
412fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4133b224e63SBarry Smith       ilink = ilink->next;
4144aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4154aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
416fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4177d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
418519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
419a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4201cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
421e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
422a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
423a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
4243b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4253b224e63SBarry Smith 
4263b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4279005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4281cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
429084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
430084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
431084e4875SJed Brown         PC pc;
432cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
433084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
434084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
435e69d4d44SBarry Smith       }
4366c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4376c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4383b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4393b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4403b224e63SBarry Smith 
4413b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4423b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4433b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4443b224e63SBarry Smith       ilink = jac->head;
4453b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4463b224e63SBarry Smith       ilink = ilink->next;
4473b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4483b224e63SBarry Smith     }
4493b224e63SBarry Smith   } else {
45097bbdb24SBarry Smith     /* set up the individual PCs */
45197bbdb24SBarry Smith     i    = 0;
4525a9f2f41SSatish Balay     ilink = jac->head;
4535a9f2f41SSatish Balay     while (ilink) {
454519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
4553b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4565a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
4575a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
45897bbdb24SBarry Smith       i++;
4595a9f2f41SSatish Balay       ilink = ilink->next;
4600971522cSBarry Smith     }
46197bbdb24SBarry Smith 
46297bbdb24SBarry Smith     /* create work vectors for each split */
4631b9fc7fcSBarry Smith     if (!jac->x) {
46497bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
4655a9f2f41SSatish Balay       ilink = jac->head;
46697bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
467906ed7ccSBarry Smith         Vec *vl,*vr;
468906ed7ccSBarry Smith 
469906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
470906ed7ccSBarry Smith         ilink->x  = *vr;
471906ed7ccSBarry Smith         ilink->y  = *vl;
472906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
473906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
4745a9f2f41SSatish Balay         jac->x[i] = ilink->x;
4755a9f2f41SSatish Balay         jac->y[i] = ilink->y;
4765a9f2f41SSatish Balay         ilink     = ilink->next;
47797bbdb24SBarry Smith       }
4783b224e63SBarry Smith     }
4793b224e63SBarry Smith   }
4803b224e63SBarry Smith 
4813b224e63SBarry Smith 
482d0f46423SBarry Smith   if (!jac->head->sctx) {
4833b224e63SBarry Smith     Vec xtmp;
4843b224e63SBarry Smith 
48579416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
4861b9fc7fcSBarry Smith 
4875a9f2f41SSatish Balay     ilink = jac->head;
4881b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
4891b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
490704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
4915a9f2f41SSatish Balay       ilink = ilink->next;
4921b9fc7fcSBarry Smith     }
4936bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
4941b9fc7fcSBarry Smith   }
4950971522cSBarry Smith   PetscFunctionReturn(0);
4960971522cSBarry Smith }
4970971522cSBarry Smith 
4985a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
499ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
500ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5015a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
502ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
503ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
50479416396SBarry Smith 
5050971522cSBarry Smith #undef __FUNCT__
5063b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5073b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5083b224e63SBarry Smith {
5093b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5103b224e63SBarry Smith   PetscErrorCode    ierr;
5113b224e63SBarry Smith   KSP               ksp;
5123b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5133b224e63SBarry Smith 
5143b224e63SBarry Smith   PetscFunctionBegin;
5153b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5163b224e63SBarry Smith 
517c5d2311dSJed Brown   switch (jac->schurfactorization) {
518c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
519a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
520c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
521c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
522c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
523c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);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,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
526c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
527c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
528c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
529c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
531c5d2311dSJed Brown       break;
532c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
533a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
534c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
535c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
536c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
537c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
538c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
539c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
540c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
541c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
542c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
543c5d2311dSJed Brown       ierr = VecScatterBegin(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       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
546c5d2311dSJed Brown       break;
547c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
548a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
549c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
550c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
551c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
552c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
553c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
554c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
555c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
556c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
557c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
558c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
559c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
560c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
561c5d2311dSJed Brown       break;
562c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
563a04f6461SBarry 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 */
5643b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5653b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5663b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5673b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
5683b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
5693b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5703b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5713b224e63SBarry Smith 
5723b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
5733b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5743b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5753b224e63SBarry Smith 
5763b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
5773b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
5783b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5793b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5803b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
581c5d2311dSJed Brown   }
5823b224e63SBarry Smith   PetscFunctionReturn(0);
5833b224e63SBarry Smith }
5843b224e63SBarry Smith 
5853b224e63SBarry Smith #undef __FUNCT__
5860971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
5870971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
5880971522cSBarry Smith {
5890971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5900971522cSBarry Smith   PetscErrorCode    ierr;
5915a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
592af93645bSJed Brown   PetscInt          cnt;
5930971522cSBarry Smith 
5940971522cSBarry Smith   PetscFunctionBegin;
59551f519a2SBarry Smith   CHKMEMQ;
59651f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
59751f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
59851f519a2SBarry Smith 
59979416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6001b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6010971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6025a9f2f41SSatish Balay       while (ilink) {
6035a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6045a9f2f41SSatish Balay         ilink = ilink->next;
6050971522cSBarry Smith       }
6060971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6071b9fc7fcSBarry Smith     } else {
608efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6095a9f2f41SSatish Balay       while (ilink) {
6105a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6115a9f2f41SSatish Balay         ilink = ilink->next;
6121b9fc7fcSBarry Smith       }
6131b9fc7fcSBarry Smith     }
61416913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
61579416396SBarry Smith     if (!jac->w1) {
61679416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
61779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
61879416396SBarry Smith     }
619efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6205a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6213e197d65SBarry Smith     cnt = 1;
6225a9f2f41SSatish Balay     while (ilink->next) {
6235a9f2f41SSatish Balay       ilink = ilink->next;
6243e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6253e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6263e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6273e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6283e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6293e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6303e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6313e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6323e197d65SBarry Smith     }
63351f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
63411755939SBarry Smith       cnt -= 2;
63551f519a2SBarry Smith       while (ilink->previous) {
63651f519a2SBarry Smith         ilink = ilink->previous;
63711755939SBarry Smith         /* compute the residual only over the part of the vector needed */
63811755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
63911755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
64011755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64111755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64211755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
64311755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64411755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64551f519a2SBarry Smith       }
64611755939SBarry Smith     }
64765e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
64851f519a2SBarry Smith   CHKMEMQ;
6490971522cSBarry Smith   PetscFunctionReturn(0);
6500971522cSBarry Smith }
6510971522cSBarry Smith 
652421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
653ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
654ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
655421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
656ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
657ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
658421e10b8SBarry Smith 
659421e10b8SBarry Smith #undef __FUNCT__
6608c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
661421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
662421e10b8SBarry Smith {
663421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
664421e10b8SBarry Smith   PetscErrorCode    ierr;
665421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
666421e10b8SBarry Smith 
667421e10b8SBarry Smith   PetscFunctionBegin;
668421e10b8SBarry Smith   CHKMEMQ;
669421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
670421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
671421e10b8SBarry Smith 
672421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
673421e10b8SBarry Smith     if (jac->defaultsplit) {
674421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
675421e10b8SBarry Smith       while (ilink) {
676421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
677421e10b8SBarry Smith 	ilink = ilink->next;
678421e10b8SBarry Smith       }
679421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
680421e10b8SBarry Smith     } else {
681421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
682421e10b8SBarry Smith       while (ilink) {
683421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
684421e10b8SBarry Smith 	ilink = ilink->next;
685421e10b8SBarry Smith       }
686421e10b8SBarry Smith     }
687421e10b8SBarry Smith   } else {
688421e10b8SBarry Smith     if (!jac->w1) {
689421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
690421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
691421e10b8SBarry Smith     }
692421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
693421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
694421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
695421e10b8SBarry Smith       while (ilink->next) {
696421e10b8SBarry Smith         ilink = ilink->next;
6979989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
698421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
699421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
700421e10b8SBarry Smith       }
701421e10b8SBarry Smith       while (ilink->previous) {
702421e10b8SBarry Smith         ilink = ilink->previous;
7039989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
704421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
705421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
706421e10b8SBarry Smith       }
707421e10b8SBarry Smith     } else {
708421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
709421e10b8SBarry Smith 	ilink = ilink->next;
710421e10b8SBarry Smith       }
711421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
712421e10b8SBarry Smith       while (ilink->previous) {
713421e10b8SBarry Smith 	ilink = ilink->previous;
7149989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
715421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
716421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
717421e10b8SBarry Smith       }
718421e10b8SBarry Smith     }
719421e10b8SBarry Smith   }
720421e10b8SBarry Smith   CHKMEMQ;
721421e10b8SBarry Smith   PetscFunctionReturn(0);
722421e10b8SBarry Smith }
723421e10b8SBarry Smith 
7240971522cSBarry Smith #undef __FUNCT__
725574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
726574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
7270971522cSBarry Smith {
7280971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7290971522cSBarry Smith   PetscErrorCode    ierr;
7305a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7310971522cSBarry Smith 
7320971522cSBarry Smith   PetscFunctionBegin;
7335a9f2f41SSatish Balay   while (ilink) {
734574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
735fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
736fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
737fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
738fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
7395a9f2f41SSatish Balay     next = ilink->next;
7405a9f2f41SSatish Balay     ilink = next;
7410971522cSBarry Smith   }
74205b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
743574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
744574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
745574deadeSBarry Smith   } else if (jac->mat) {
746574deadeSBarry Smith     jac->mat = PETSC_NULL;
747574deadeSBarry Smith   }
74897bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
74968dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
7506bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
7516bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
7526bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
7536bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
7546bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
7556bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
7566bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
757574deadeSBarry Smith   PetscFunctionReturn(0);
758574deadeSBarry Smith }
759574deadeSBarry Smith 
760574deadeSBarry Smith #undef __FUNCT__
761574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
762574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
763574deadeSBarry Smith {
764574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
765574deadeSBarry Smith   PetscErrorCode    ierr;
766574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
767574deadeSBarry Smith 
768574deadeSBarry Smith   PetscFunctionBegin;
769574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
770574deadeSBarry Smith   while (ilink) {
7716bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
772574deadeSBarry Smith     next = ilink->next;
773574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
774574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
775574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
776574deadeSBarry Smith     ilink = next;
777574deadeSBarry Smith   }
778574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
779c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
7800971522cSBarry Smith   PetscFunctionReturn(0);
7810971522cSBarry Smith }
7820971522cSBarry Smith 
7830971522cSBarry Smith #undef __FUNCT__
7840971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
7850971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
7860971522cSBarry Smith {
7871b9fc7fcSBarry Smith   PetscErrorCode  ierr;
7886c924f48SJed Brown   PetscInt        bs;
789bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
7909dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
7913b224e63SBarry Smith   PCCompositeType ctype;
7921b9fc7fcSBarry Smith 
7930971522cSBarry Smith   PetscFunctionBegin;
7941b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
795acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
79651f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
79751f519a2SBarry Smith   if (flg) {
79851f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
79951f519a2SBarry Smith   }
800704ba839SBarry Smith 
801435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
802c0adfefeSBarry Smith   if (stokes) {
803c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
804c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
805c0adfefeSBarry Smith   }
806c0adfefeSBarry Smith 
8073b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8083b224e63SBarry Smith   if (flg) {
8093b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8103b224e63SBarry Smith   }
811d32f9abdSBarry Smith 
812c30613efSMatthew Knepley   /* Only setup fields once */
813c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
814d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
815d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8166c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8176c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
818d32f9abdSBarry Smith   }
819c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
820c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
821084e4875SJed 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);
822c5d2311dSJed Brown   }
8231b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8240971522cSBarry Smith   PetscFunctionReturn(0);
8250971522cSBarry Smith }
8260971522cSBarry Smith 
8270971522cSBarry Smith /*------------------------------------------------------------------------------------*/
8280971522cSBarry Smith 
8290971522cSBarry Smith EXTERN_C_BEGIN
8300971522cSBarry Smith #undef __FUNCT__
8310971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
8327087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
8330971522cSBarry Smith {
83497bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8350971522cSBarry Smith   PetscErrorCode    ierr;
8365a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
83769a612a9SBarry Smith   char              prefix[128];
83851f519a2SBarry Smith   PetscInt          i;
8390971522cSBarry Smith 
8400971522cSBarry Smith   PetscFunctionBegin;
8416c924f48SJed Brown   if (jac->splitdefined) {
8426c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8436c924f48SJed Brown     PetscFunctionReturn(0);
8446c924f48SJed Brown   }
84551f519a2SBarry Smith   for (i=0; i<n; i++) {
846e32f2f54SBarry 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);
847e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
84851f519a2SBarry Smith   }
849704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
850a04f6461SBarry Smith   if (splitname) {
851db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
852a04f6461SBarry Smith   } else {
853a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
854a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
855a04f6461SBarry Smith   }
856704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
8575a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
8585a9f2f41SSatish Balay   ilink->nfields = n;
8595a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
8607adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
8611cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
8625a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
8639005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
86469a612a9SBarry Smith 
865a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
8665a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
8670971522cSBarry Smith 
8680971522cSBarry Smith   if (!next) {
8695a9f2f41SSatish Balay     jac->head       = ilink;
87051f519a2SBarry Smith     ilink->previous = PETSC_NULL;
8710971522cSBarry Smith   } else {
8720971522cSBarry Smith     while (next->next) {
8730971522cSBarry Smith       next = next->next;
8740971522cSBarry Smith     }
8755a9f2f41SSatish Balay     next->next      = ilink;
87651f519a2SBarry Smith     ilink->previous = next;
8770971522cSBarry Smith   }
8780971522cSBarry Smith   jac->nsplits++;
8790971522cSBarry Smith   PetscFunctionReturn(0);
8800971522cSBarry Smith }
8810971522cSBarry Smith EXTERN_C_END
8820971522cSBarry Smith 
883e69d4d44SBarry Smith EXTERN_C_BEGIN
884e69d4d44SBarry Smith #undef __FUNCT__
885e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
8867087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
887e69d4d44SBarry Smith {
888e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
889e69d4d44SBarry Smith   PetscErrorCode ierr;
890e69d4d44SBarry Smith 
891e69d4d44SBarry Smith   PetscFunctionBegin;
892519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
893e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
894e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
89513e0d083SBarry Smith   if (n) *n = jac->nsplits;
896e69d4d44SBarry Smith   PetscFunctionReturn(0);
897e69d4d44SBarry Smith }
898e69d4d44SBarry Smith EXTERN_C_END
8990971522cSBarry Smith 
9000971522cSBarry Smith EXTERN_C_BEGIN
9010971522cSBarry Smith #undef __FUNCT__
90269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9037087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9040971522cSBarry Smith {
9050971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9060971522cSBarry Smith   PetscErrorCode    ierr;
9070971522cSBarry Smith   PetscInt          cnt = 0;
9085a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9090971522cSBarry Smith 
9100971522cSBarry Smith   PetscFunctionBegin;
911*5d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9125a9f2f41SSatish Balay   while (ilink) {
9135a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9145a9f2f41SSatish Balay     ilink = ilink->next;
9150971522cSBarry Smith   }
916*5d480477SMatthew G Knepley   if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number of splits in linked list %D does not match number in object %D",cnt,jac->nsplits);
91713e0d083SBarry Smith   if (n) *n = jac->nsplits;
9180971522cSBarry Smith   PetscFunctionReturn(0);
9190971522cSBarry Smith }
9200971522cSBarry Smith EXTERN_C_END
9210971522cSBarry Smith 
922704ba839SBarry Smith EXTERN_C_BEGIN
923704ba839SBarry Smith #undef __FUNCT__
924704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
9257087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
926704ba839SBarry Smith {
927704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
928704ba839SBarry Smith   PetscErrorCode    ierr;
929704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
930704ba839SBarry Smith   char              prefix[128];
931704ba839SBarry Smith 
932704ba839SBarry Smith   PetscFunctionBegin;
9336c924f48SJed Brown   if (jac->splitdefined) {
9346c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9356c924f48SJed Brown     PetscFunctionReturn(0);
9366c924f48SJed Brown   }
93716913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
938a04f6461SBarry Smith   if (splitname) {
939db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
940a04f6461SBarry Smith   } else {
941a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
942a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
943a04f6461SBarry Smith   }
9441ab39975SBarry Smith   ilink->is      = is;
9451ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
946704ba839SBarry Smith   ilink->next    = PETSC_NULL;
947704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9481cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
949704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9509005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
951704ba839SBarry Smith 
952a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
953704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
954704ba839SBarry Smith 
955704ba839SBarry Smith   if (!next) {
956704ba839SBarry Smith     jac->head       = ilink;
957704ba839SBarry Smith     ilink->previous = PETSC_NULL;
958704ba839SBarry Smith   } else {
959704ba839SBarry Smith     while (next->next) {
960704ba839SBarry Smith       next = next->next;
961704ba839SBarry Smith     }
962704ba839SBarry Smith     next->next      = ilink;
963704ba839SBarry Smith     ilink->previous = next;
964704ba839SBarry Smith   }
965704ba839SBarry Smith   jac->nsplits++;
966704ba839SBarry Smith 
967704ba839SBarry Smith   PetscFunctionReturn(0);
968704ba839SBarry Smith }
969704ba839SBarry Smith EXTERN_C_END
970704ba839SBarry Smith 
9710971522cSBarry Smith #undef __FUNCT__
9720971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
9730971522cSBarry Smith /*@
9740971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
9750971522cSBarry Smith 
976ad4df100SBarry Smith     Logically Collective on PC
9770971522cSBarry Smith 
9780971522cSBarry Smith     Input Parameters:
9790971522cSBarry Smith +   pc  - the preconditioner context
980a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
9810971522cSBarry Smith .   n - the number of fields in this split
982db4c96c1SJed Brown -   fields - the fields in this split
9830971522cSBarry Smith 
9840971522cSBarry Smith     Level: intermediate
9850971522cSBarry Smith 
986d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
987d32f9abdSBarry Smith 
988d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
989d32f9abdSBarry 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
990d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
991d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
992d32f9abdSBarry Smith 
993db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
994db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
995db4c96c1SJed Brown 
996d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
9970971522cSBarry Smith 
9980971522cSBarry Smith @*/
9997087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10000971522cSBarry Smith {
10014ac538c5SBarry Smith   PetscErrorCode ierr;
10020971522cSBarry Smith 
10030971522cSBarry Smith   PetscFunctionBegin;
10040700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1005db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1006db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1007db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10084ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10090971522cSBarry Smith   PetscFunctionReturn(0);
10100971522cSBarry Smith }
10110971522cSBarry Smith 
10120971522cSBarry Smith #undef __FUNCT__
1013704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1014704ba839SBarry Smith /*@
1015704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1016704ba839SBarry Smith 
1017ad4df100SBarry Smith     Logically Collective on PC
1018704ba839SBarry Smith 
1019704ba839SBarry Smith     Input Parameters:
1020704ba839SBarry Smith +   pc  - the preconditioner context
1021a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1022db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1023704ba839SBarry Smith 
1024d32f9abdSBarry Smith 
1025a6ffb8dbSJed Brown     Notes:
1026a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1027a6ffb8dbSJed Brown 
1028db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1029db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1030d32f9abdSBarry Smith 
1031704ba839SBarry Smith     Level: intermediate
1032704ba839SBarry Smith 
1033704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1034704ba839SBarry Smith 
1035704ba839SBarry Smith @*/
10367087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1037704ba839SBarry Smith {
10384ac538c5SBarry Smith   PetscErrorCode ierr;
1039704ba839SBarry Smith 
1040704ba839SBarry Smith   PetscFunctionBegin;
10410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1042db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1043db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
10444ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1045704ba839SBarry Smith   PetscFunctionReturn(0);
1046704ba839SBarry Smith }
1047704ba839SBarry Smith 
1048704ba839SBarry Smith #undef __FUNCT__
104957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
105057a9adfeSMatthew G Knepley /*@
105157a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
105257a9adfeSMatthew G Knepley 
105357a9adfeSMatthew G Knepley     Logically Collective on PC
105457a9adfeSMatthew G Knepley 
105557a9adfeSMatthew G Knepley     Input Parameters:
105657a9adfeSMatthew G Knepley +   pc  - the preconditioner context
105757a9adfeSMatthew G Knepley -   splitname - name of this split
105857a9adfeSMatthew G Knepley 
105957a9adfeSMatthew G Knepley     Output Parameter:
106057a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
106157a9adfeSMatthew G Knepley 
106257a9adfeSMatthew G Knepley     Level: intermediate
106357a9adfeSMatthew G Knepley 
106457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
106557a9adfeSMatthew G Knepley 
106657a9adfeSMatthew G Knepley @*/
106757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
106857a9adfeSMatthew G Knepley {
106957a9adfeSMatthew G Knepley   PetscErrorCode ierr;
107057a9adfeSMatthew G Knepley 
107157a9adfeSMatthew G Knepley   PetscFunctionBegin;
107257a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
107357a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
107457a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
107557a9adfeSMatthew G Knepley   {
107657a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
107757a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
107857a9adfeSMatthew G Knepley     PetscBool         found;
107957a9adfeSMatthew G Knepley 
108057a9adfeSMatthew G Knepley     *is = PETSC_NULL;
108157a9adfeSMatthew G Knepley     while(ilink) {
108257a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
108357a9adfeSMatthew G Knepley       if (found) {
108457a9adfeSMatthew G Knepley         *is = ilink->is;
108557a9adfeSMatthew G Knepley         break;
108657a9adfeSMatthew G Knepley       }
108757a9adfeSMatthew G Knepley       ilink = ilink->next;
108857a9adfeSMatthew G Knepley     }
108957a9adfeSMatthew G Knepley   }
109057a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
109157a9adfeSMatthew G Knepley }
109257a9adfeSMatthew G Knepley 
109357a9adfeSMatthew G Knepley #undef __FUNCT__
109451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
109551f519a2SBarry Smith /*@
109651f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
109751f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
109851f519a2SBarry Smith 
1099ad4df100SBarry Smith     Logically Collective on PC
110051f519a2SBarry Smith 
110151f519a2SBarry Smith     Input Parameters:
110251f519a2SBarry Smith +   pc  - the preconditioner context
110351f519a2SBarry Smith -   bs - the block size
110451f519a2SBarry Smith 
110551f519a2SBarry Smith     Level: intermediate
110651f519a2SBarry Smith 
110751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
110851f519a2SBarry Smith 
110951f519a2SBarry Smith @*/
11107087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
111151f519a2SBarry Smith {
11124ac538c5SBarry Smith   PetscErrorCode ierr;
111351f519a2SBarry Smith 
111451f519a2SBarry Smith   PetscFunctionBegin;
11150700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1116c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
11174ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
111851f519a2SBarry Smith   PetscFunctionReturn(0);
111951f519a2SBarry Smith }
112051f519a2SBarry Smith 
112151f519a2SBarry Smith #undef __FUNCT__
112269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
11230971522cSBarry Smith /*@C
112469a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
11250971522cSBarry Smith 
112669a612a9SBarry Smith    Collective on KSP
11270971522cSBarry Smith 
11280971522cSBarry Smith    Input Parameter:
11290971522cSBarry Smith .  pc - the preconditioner context
11300971522cSBarry Smith 
11310971522cSBarry Smith    Output Parameters:
113213e0d083SBarry Smith +  n - the number of splits
113369a612a9SBarry Smith -  pc - the array of KSP contexts
11340971522cSBarry Smith 
11350971522cSBarry Smith    Note:
1136d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1137d32f9abdSBarry Smith    (not the KSP just the array that contains them).
11380971522cSBarry Smith 
113969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
11400971522cSBarry Smith 
11410971522cSBarry Smith    Level: advanced
11420971522cSBarry Smith 
11430971522cSBarry Smith .seealso: PCFIELDSPLIT
11440971522cSBarry Smith @*/
11457087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
11460971522cSBarry Smith {
11474ac538c5SBarry Smith   PetscErrorCode ierr;
11480971522cSBarry Smith 
11490971522cSBarry Smith   PetscFunctionBegin;
11500700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115113e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
11524ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
11530971522cSBarry Smith   PetscFunctionReturn(0);
11540971522cSBarry Smith }
11550971522cSBarry Smith 
1156e69d4d44SBarry Smith #undef __FUNCT__
1157e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1158e69d4d44SBarry Smith /*@
1159e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1160a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1161e69d4d44SBarry Smith 
1162e69d4d44SBarry Smith     Collective on PC
1163e69d4d44SBarry Smith 
1164e69d4d44SBarry Smith     Input Parameters:
1165e69d4d44SBarry Smith +   pc  - the preconditioner context
1166084e4875SJed Brown .   ptype - which matrix to use for preconditioning the Schur complement
1167084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1168084e4875SJed Brown 
1169084e4875SJed Brown     Notes:
1170a04f6461SBarry Smith     The default is to use the block on the diagonal of the preconditioning matrix.  This is A11, in the (1,1) position.
1171084e4875SJed Brown     There are currently no preconditioners that work directly with the Schur complement so setting
1172084e4875SJed Brown     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
1173e69d4d44SBarry Smith 
1174e69d4d44SBarry Smith     Options Database:
1175084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1176e69d4d44SBarry Smith 
1177e69d4d44SBarry Smith     Level: intermediate
1178e69d4d44SBarry Smith 
1179084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1180e69d4d44SBarry Smith 
1181e69d4d44SBarry Smith @*/
11827087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1183e69d4d44SBarry Smith {
11844ac538c5SBarry Smith   PetscErrorCode ierr;
1185e69d4d44SBarry Smith 
1186e69d4d44SBarry Smith   PetscFunctionBegin;
11870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11884ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1189e69d4d44SBarry Smith   PetscFunctionReturn(0);
1190e69d4d44SBarry Smith }
1191e69d4d44SBarry Smith 
1192e69d4d44SBarry Smith EXTERN_C_BEGIN
1193e69d4d44SBarry Smith #undef __FUNCT__
1194e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
11957087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1196e69d4d44SBarry Smith {
1197e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1198084e4875SJed Brown   PetscErrorCode  ierr;
1199e69d4d44SBarry Smith 
1200e69d4d44SBarry Smith   PetscFunctionBegin;
1201084e4875SJed Brown   jac->schurpre = ptype;
1202084e4875SJed Brown   if (pre) {
12036bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1204084e4875SJed Brown     jac->schur_user = pre;
1205084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1206084e4875SJed Brown   }
1207e69d4d44SBarry Smith   PetscFunctionReturn(0);
1208e69d4d44SBarry Smith }
1209e69d4d44SBarry Smith EXTERN_C_END
1210e69d4d44SBarry Smith 
121130ad9308SMatthew Knepley #undef __FUNCT__
121230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
121330ad9308SMatthew Knepley /*@C
12148c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
121530ad9308SMatthew Knepley 
121630ad9308SMatthew Knepley    Collective on KSP
121730ad9308SMatthew Knepley 
121830ad9308SMatthew Knepley    Input Parameter:
121930ad9308SMatthew Knepley .  pc - the preconditioner context
122030ad9308SMatthew Knepley 
122130ad9308SMatthew Knepley    Output Parameters:
1222a04f6461SBarry Smith +  A00 - the (0,0) block
1223a04f6461SBarry Smith .  A01 - the (0,1) block
1224a04f6461SBarry Smith .  A10 - the (1,0) block
1225a04f6461SBarry Smith -  A11 - the (1,1) block
122630ad9308SMatthew Knepley 
122730ad9308SMatthew Knepley    Level: advanced
122830ad9308SMatthew Knepley 
122930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
123030ad9308SMatthew Knepley @*/
1231a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
123230ad9308SMatthew Knepley {
123330ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
123430ad9308SMatthew Knepley 
123530ad9308SMatthew Knepley   PetscFunctionBegin;
12360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1237c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1238a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1239a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1240a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1241a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
124230ad9308SMatthew Knepley   PetscFunctionReturn(0);
124330ad9308SMatthew Knepley }
124430ad9308SMatthew Knepley 
124579416396SBarry Smith EXTERN_C_BEGIN
124679416396SBarry Smith #undef __FUNCT__
124779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
12487087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
124979416396SBarry Smith {
125079416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1251e69d4d44SBarry Smith   PetscErrorCode ierr;
125279416396SBarry Smith 
125379416396SBarry Smith   PetscFunctionBegin;
125479416396SBarry Smith   jac->type = type;
12553b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
12563b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
12573b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1258e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1259e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1260e69d4d44SBarry Smith 
12613b224e63SBarry Smith   } else {
12623b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
12633b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1264e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
12659e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
12663b224e63SBarry Smith   }
126779416396SBarry Smith   PetscFunctionReturn(0);
126879416396SBarry Smith }
126979416396SBarry Smith EXTERN_C_END
127079416396SBarry Smith 
127151f519a2SBarry Smith EXTERN_C_BEGIN
127251f519a2SBarry Smith #undef __FUNCT__
127351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
12747087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
127551f519a2SBarry Smith {
127651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
127751f519a2SBarry Smith 
127851f519a2SBarry Smith   PetscFunctionBegin;
127965e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
128065e19b50SBarry 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);
128151f519a2SBarry Smith   jac->bs = bs;
128251f519a2SBarry Smith   PetscFunctionReturn(0);
128351f519a2SBarry Smith }
128451f519a2SBarry Smith EXTERN_C_END
128551f519a2SBarry Smith 
128679416396SBarry Smith #undef __FUNCT__
128779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1288bc08b0f1SBarry Smith /*@
128979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
129079416396SBarry Smith 
129179416396SBarry Smith    Collective on PC
129279416396SBarry Smith 
129379416396SBarry Smith    Input Parameter:
129479416396SBarry Smith .  pc - the preconditioner context
129581540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
129679416396SBarry Smith 
129779416396SBarry Smith    Options Database Key:
1298a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
129979416396SBarry Smith 
130079416396SBarry Smith    Level: Developer
130179416396SBarry Smith 
130279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
130379416396SBarry Smith 
130479416396SBarry Smith .seealso: PCCompositeSetType()
130579416396SBarry Smith 
130679416396SBarry Smith @*/
13077087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
130879416396SBarry Smith {
13094ac538c5SBarry Smith   PetscErrorCode ierr;
131079416396SBarry Smith 
131179416396SBarry Smith   PetscFunctionBegin;
13120700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13134ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
131479416396SBarry Smith   PetscFunctionReturn(0);
131579416396SBarry Smith }
131679416396SBarry Smith 
13170971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
13180971522cSBarry Smith /*MC
1319a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1320a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
13210971522cSBarry Smith 
1322edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1323edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
13240971522cSBarry Smith 
1325a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
132669a612a9SBarry Smith          and set the options directly on the resulting KSP object
13270971522cSBarry Smith 
13280971522cSBarry Smith    Level: intermediate
13290971522cSBarry Smith 
133079416396SBarry Smith    Options Database Keys:
133181540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
133281540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
133381540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
133481540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
133581540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1336e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1337435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
133879416396SBarry Smith 
1339edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
13403b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
13413b224e63SBarry Smith 
1342d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1343d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1344d32f9abdSBarry Smith 
1345d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1346d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1347d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1348d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1349d32f9abdSBarry Smith 
1350a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1351a04f6461SBarry Smith                                                     ( A10 A11 )
1352e69d4d44SBarry Smith      the preconditioner is
1353a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1354a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1355a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1356a04f6461SBarry 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).
1357a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1358edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1359a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
13607e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1361e69d4d44SBarry Smith 
1362edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1363edf189efSBarry Smith      is used automatically for a second block.
1364edf189efSBarry Smith 
13650716a85fSBarry Smith      The fieldsplit preconditioner cannot be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. Generally it should be used with the AIJ format.
13660716a85fSBarry Smith 
1367a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
13680971522cSBarry Smith 
13697e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1370e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
13710971522cSBarry Smith M*/
13720971522cSBarry Smith 
13730971522cSBarry Smith EXTERN_C_BEGIN
13740971522cSBarry Smith #undef __FUNCT__
13750971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
13767087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
13770971522cSBarry Smith {
13780971522cSBarry Smith   PetscErrorCode ierr;
13790971522cSBarry Smith   PC_FieldSplit  *jac;
13800971522cSBarry Smith 
13810971522cSBarry Smith   PetscFunctionBegin;
138238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
13830971522cSBarry Smith   jac->bs        = -1;
13840971522cSBarry Smith   jac->nsplits   = 0;
13853e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1386e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1387c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
138851f519a2SBarry Smith 
13890971522cSBarry Smith   pc->data     = (void*)jac;
13900971522cSBarry Smith 
13910971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1392421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
13930971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1394574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
13950971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
13960971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
13970971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
13980971522cSBarry Smith   pc->ops->applyrichardson   = 0;
13990971522cSBarry Smith 
140069a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
140169a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14020971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
14030971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1404704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1405704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
140679416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
140779416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
140851f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
140951f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
14100971522cSBarry Smith   PetscFunctionReturn(0);
14110971522cSBarry Smith }
14120971522cSBarry Smith EXTERN_C_END
14130971522cSBarry Smith 
14140971522cSBarry Smith 
1415a541d17aSBarry Smith 
1416