xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 7ae8954a65a71528aba8abf0b6fd4109ea6196c4)
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;
4863ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
49c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
500971522cSBarry Smith } PC_FieldSplit;
510971522cSBarry Smith 
5216913363SBarry Smith /*
5316913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5416913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5516913363SBarry Smith    PC you could change this.
5616913363SBarry Smith */
57084e4875SJed Brown 
58e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
59084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
60084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
61084e4875SJed Brown {
62084e4875SJed Brown   switch (jac->schurpre) {
63084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
64084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
65084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
66084e4875SJed Brown     default:
67084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
68084e4875SJed Brown   }
69084e4875SJed Brown }
70084e4875SJed Brown 
71084e4875SJed Brown 
720971522cSBarry Smith #undef __FUNCT__
730971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
740971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
750971522cSBarry Smith {
760971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
770971522cSBarry Smith   PetscErrorCode    ierr;
78ace3abfcSBarry Smith   PetscBool         iascii;
790971522cSBarry Smith   PetscInt          i,j;
805a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
810971522cSBarry Smith 
820971522cSBarry Smith   PetscFunctionBegin;
832692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
840971522cSBarry Smith   if (iascii) {
8551f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
870971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
880971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
891ab39975SBarry Smith       if (ilink->fields) {
900971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9179416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
925a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9379416396SBarry Smith 	  if (j > 0) {
9479416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9579416396SBarry Smith 	  }
965a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
970971522cSBarry Smith 	}
980971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1001ab39975SBarry Smith       } else {
1011ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1021ab39975SBarry Smith       }
1035a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1045a9f2f41SSatish Balay       ilink = ilink->next;
1050971522cSBarry Smith     }
1060971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1070971522cSBarry Smith   } else {
10865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1090971522cSBarry Smith   }
1100971522cSBarry Smith   PetscFunctionReturn(0);
1110971522cSBarry Smith }
1120971522cSBarry Smith 
1130971522cSBarry Smith #undef __FUNCT__
1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1163b224e63SBarry Smith {
1173b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1183b224e63SBarry Smith   PetscErrorCode    ierr;
119ace3abfcSBarry Smith   PetscBool         iascii;
1203b224e63SBarry Smith   PetscInt          i,j;
1213b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1223b224e63SBarry Smith   KSP               ksp;
1233b224e63SBarry Smith 
1243b224e63SBarry Smith   PetscFunctionBegin;
1252692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1263b224e63SBarry Smith   if (iascii) {
127c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1283b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1303b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1313b224e63SBarry Smith       if (ilink->fields) {
1323b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1343b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1353b224e63SBarry Smith 	  if (j > 0) {
1363b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1373b224e63SBarry Smith 	  }
1383b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1393b224e63SBarry Smith 	}
1403b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1413b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1423b224e63SBarry Smith       } else {
1433b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1443b224e63SBarry Smith       }
1453b224e63SBarry Smith       ilink = ilink->next;
1463b224e63SBarry Smith     }
147435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1483b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14912cae6f2SJed Brown     if (jac->schur) {
15012cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1513b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15212cae6f2SJed Brown     } else {
15312cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15412cae6f2SJed Brown     }
1553b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
156435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1573b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15812cae6f2SJed Brown     if (jac->kspschur) {
1593b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16012cae6f2SJed Brown     } else {
16112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16212cae6f2SJed Brown     }
1633b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith   } else {
16665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1673b224e63SBarry Smith   }
1683b224e63SBarry Smith   PetscFunctionReturn(0);
1693b224e63SBarry Smith }
1703b224e63SBarry Smith 
1713b224e63SBarry Smith #undef __FUNCT__
1726c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1736c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1746c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1756c924f48SJed Brown {
1766c924f48SJed Brown   PetscErrorCode ierr;
1776c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1786c924f48SJed Brown   PetscInt       i,nfields,*ifields;
179ace3abfcSBarry Smith   PetscBool      flg;
1806c924f48SJed Brown   char           optionname[128],splitname[8];
1816c924f48SJed Brown 
1826c924f48SJed Brown   PetscFunctionBegin;
1836c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1846c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1856c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1866c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1876c924f48SJed Brown     nfields = jac->bs;
1886c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1896c924f48SJed Brown     if (!flg) break;
1906c924f48SJed Brown     if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1916c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1926c924f48SJed Brown   }
1936c924f48SJed Brown   if (i > 0) {
1946c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1956c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1966c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1976c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1986c924f48SJed Brown   }
1996c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2006c924f48SJed Brown   PetscFunctionReturn(0);
2016c924f48SJed Brown }
2026c924f48SJed Brown 
2036c924f48SJed Brown #undef __FUNCT__
20469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20569a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2060971522cSBarry Smith {
2070971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2080971522cSBarry Smith   PetscErrorCode    ierr;
2095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2106ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2116c924f48SJed Brown   PetscInt          i;
2120971522cSBarry Smith 
2130971522cSBarry Smith   PetscFunctionBegin;
214d32f9abdSBarry Smith   if (!ilink) {
215d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
216d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2178b8307b2SJed Brown       PetscBool dmcomposite;
2188b8307b2SJed Brown       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
2198b8307b2SJed Brown       if (dmcomposite) {
2208b8307b2SJed Brown         PetscInt nDM;
2218b8307b2SJed Brown         IS       *fields;
2222fa5ba8aSJed Brown         DM       *dms;
2238b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2248b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
2258b8307b2SJed Brown         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
2262fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
2272fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
2288b8307b2SJed Brown         for (i=0; i<nDM; i++) {
2292fa5ba8aSJed Brown           char buf[256];
2302fa5ba8aSJed Brown           const char *splitname;
2312fa5ba8aSJed Brown           /* Split naming precedence: object name, prefix, number */
2322fa5ba8aSJed Brown           splitname = ((PetscObject)pc->dm)->name;
2332fa5ba8aSJed Brown           if (!splitname) {
2342fa5ba8aSJed Brown             ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
2352fa5ba8aSJed Brown             if (splitname) {
2362fa5ba8aSJed Brown               size_t len;
2372fa5ba8aSJed Brown               ierr = PetscStrncpy(buf,splitname,sizeof buf);CHKERRQ(ierr);
2382fa5ba8aSJed Brown               buf[sizeof buf - 1] = 0;
2392fa5ba8aSJed Brown               ierr = PetscStrlen(buf,&len);CHKERRQ(ierr);
2402fa5ba8aSJed Brown               if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
2412fa5ba8aSJed Brown               splitname = buf;
2422fa5ba8aSJed Brown             }
2432fa5ba8aSJed Brown           }
2442fa5ba8aSJed Brown           if (!splitname) {
2452fa5ba8aSJed Brown             ierr = PetscSNPrintf(buf,sizeof buf,"%D",i);CHKERRQ(ierr);
2462fa5ba8aSJed Brown             splitname = buf;
2472fa5ba8aSJed Brown           }
2488b8307b2SJed Brown           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
249fcfd50ebSBarry Smith           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
2508b8307b2SJed Brown         }
2518b8307b2SJed Brown         ierr = PetscFree(fields);CHKERRQ(ierr);
2522fa5ba8aSJed Brown         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2532fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
2542fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
2552fa5ba8aSJed Brown         }
2562fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
2578b8307b2SJed Brown       }
25866ffff09SJed Brown     } else {
259521d7252SBarry Smith       if (jac->bs <= 0) {
260704ba839SBarry Smith         if (pc->pmat) {
261521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
262704ba839SBarry Smith         } else {
263704ba839SBarry Smith           jac->bs = 1;
264704ba839SBarry Smith         }
265521d7252SBarry Smith       }
266d32f9abdSBarry Smith 
267acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2686ce1633cSBarry Smith       if (stokes) {
2696ce1633cSBarry Smith         IS       zerodiags,rest;
2706ce1633cSBarry Smith         PetscInt nmin,nmax;
2716ce1633cSBarry Smith 
2726ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2736ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2746ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2756ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2766ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
277fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
278fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2796ce1633cSBarry Smith       } else {
280d32f9abdSBarry Smith         if (!flg) {
281d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
282d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2836c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2846c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
285d32f9abdSBarry Smith         }
2866c924f48SJed Brown         if (flg || !jac->splitdefined) {
287d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
288db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2896c924f48SJed Brown             char splitname[8];
2906c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
291db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
29279416396SBarry Smith           }
29397bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
294521d7252SBarry Smith         }
29566ffff09SJed Brown       }
2966ce1633cSBarry Smith     }
297edf189efSBarry Smith   } else if (jac->nsplits == 1) {
298edf189efSBarry Smith     if (ilink->is) {
299edf189efSBarry Smith       IS       is2;
300edf189efSBarry Smith       PetscInt nmin,nmax;
301edf189efSBarry Smith 
302edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
303edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
304db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
305fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
306e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
30763ec74ffSBarry Smith   } else if (jac->reset) {
30863ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
309d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
310d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
311d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
312d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
313d0af7cd3SBarry Smith       PetscBool dmcomposite;
314d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
315d0af7cd3SBarry Smith       if (dmcomposite) {
316d0af7cd3SBarry Smith         PetscInt nDM;
317d0af7cd3SBarry Smith         IS       *fields;
3182fa5ba8aSJed Brown         DM       *dms;
319d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
320d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
321d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3222fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3232fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
324d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3252fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3262fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
327d0af7cd3SBarry Smith           ilink->is = fields[i];
328d0af7cd3SBarry Smith           ilink     = ilink->next;
329edf189efSBarry Smith         }
330d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3312fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
332d0af7cd3SBarry Smith       }
333d0af7cd3SBarry Smith     } else {
334d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
335d0af7cd3SBarry Smith       if (stokes) {
336d0af7cd3SBarry Smith         IS       zerodiags,rest;
337d0af7cd3SBarry Smith         PetscInt nmin,nmax;
338d0af7cd3SBarry Smith 
339d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
340d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
341d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
34220b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
34320b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
344d0af7cd3SBarry Smith         ilink->is       = rest;
345d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
34620b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
347d0af7cd3SBarry Smith     }
348d0af7cd3SBarry Smith   }
349d0af7cd3SBarry Smith 
350e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
35169a612a9SBarry Smith   PetscFunctionReturn(0);
35269a612a9SBarry Smith }
35369a612a9SBarry Smith 
35469a612a9SBarry Smith #undef __FUNCT__
35569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
35669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
35769a612a9SBarry Smith {
35869a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
35969a612a9SBarry Smith   PetscErrorCode    ierr;
3605a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
361cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
36269a612a9SBarry Smith   MatStructure      flag = pc->flag;
363ace3abfcSBarry Smith   PetscBool         sorted;
36469a612a9SBarry Smith 
36569a612a9SBarry Smith   PetscFunctionBegin;
36669a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
36797bbdb24SBarry Smith   nsplit = jac->nsplits;
3685a9f2f41SSatish Balay   ilink  = jac->head;
36997bbdb24SBarry Smith 
37097bbdb24SBarry Smith   /* get the matrices for each split */
371704ba839SBarry Smith   if (!jac->issetup) {
3721b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
37397bbdb24SBarry Smith 
374704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
375704ba839SBarry Smith 
376704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
37751f519a2SBarry Smith     bs     = jac->bs;
37897bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
379cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3801b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3811b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3821b9fc7fcSBarry Smith       if (jac->defaultsplit) {
383704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
384704ba839SBarry Smith       } else if (!ilink->is) {
385ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3865a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3875a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3881b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3891b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3901b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
39197bbdb24SBarry Smith             }
39297bbdb24SBarry Smith           }
393d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
394ccb205f8SBarry Smith         } else {
395704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
396ccb205f8SBarry Smith         }
3973e197d65SBarry Smith       }
398704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
399e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
400704ba839SBarry Smith       ilink = ilink->next;
4011b9fc7fcSBarry Smith     }
4021b9fc7fcSBarry Smith   }
4031b9fc7fcSBarry Smith 
404704ba839SBarry Smith   ilink  = jac->head;
40597bbdb24SBarry Smith   if (!jac->pmat) {
406cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
407cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4084aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
409704ba839SBarry Smith       ilink = ilink->next;
410cf502942SBarry Smith     }
41197bbdb24SBarry Smith   } else {
412cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4134aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
414704ba839SBarry Smith       ilink = ilink->next;
415cf502942SBarry Smith     }
41697bbdb24SBarry Smith   }
417519d70e2SJed Brown   if (jac->realdiagonal) {
418519d70e2SJed Brown     ilink = jac->head;
419519d70e2SJed Brown     if (!jac->mat) {
420519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
421519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
422519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
423519d70e2SJed Brown         ilink = ilink->next;
424519d70e2SJed Brown       }
425519d70e2SJed Brown     } else {
426519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
427966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
428519d70e2SJed Brown         ilink = ilink->next;
429519d70e2SJed Brown       }
430519d70e2SJed Brown     }
431519d70e2SJed Brown   } else {
432519d70e2SJed Brown     jac->mat = jac->pmat;
433519d70e2SJed Brown   }
43497bbdb24SBarry Smith 
4356c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
43668dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
43768dd23aaSBarry Smith     ilink  = jac->head;
43868dd23aaSBarry Smith     if (!jac->Afield) {
43968dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
44068dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4414aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
44268dd23aaSBarry Smith         ilink = ilink->next;
44368dd23aaSBarry Smith       }
44468dd23aaSBarry Smith     } else {
44568dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4464aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
44768dd23aaSBarry Smith         ilink = ilink->next;
44868dd23aaSBarry Smith       }
44968dd23aaSBarry Smith     }
45068dd23aaSBarry Smith   }
45168dd23aaSBarry Smith 
4523b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4533b224e63SBarry Smith     IS       ccis;
4544aa3045dSJed Brown     PetscInt rstart,rend;
455e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
45668dd23aaSBarry Smith 
457e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
458e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
459e6cab6aaSJed Brown 
4603b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4613b224e63SBarry Smith     if (jac->schur) {
4623b224e63SBarry Smith       ilink = jac->head;
4634aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4644aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
465fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4663b224e63SBarry Smith       ilink = ilink->next;
4674aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4684aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
469fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
470519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
471084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4723b224e63SBarry Smith 
4733b224e63SBarry Smith      } else {
4741cee3971SBarry Smith       KSP ksp;
4756c924f48SJed Brown       char schurprefix[256];
4763b224e63SBarry Smith 
477a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4783b224e63SBarry Smith       ilink = jac->head;
4794aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4804aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
481fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4823b224e63SBarry Smith       ilink = ilink->next;
4834aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4844aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
485fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4867d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
487519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
488a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4891cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
490e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
491a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
492a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
49320b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
49420b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
495*7ae8954aSSatish Balay       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
4963b224e63SBarry Smith 
4973b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4989005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4991cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
500084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
501084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
502084e4875SJed Brown         PC pc;
503cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
504084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
505084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
506e69d4d44SBarry Smith       }
5076c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5086c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5093b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
51020b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
51120b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5123b224e63SBarry Smith 
5133b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5143b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5153b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5163b224e63SBarry Smith       ilink = jac->head;
5173b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5183b224e63SBarry Smith       ilink = ilink->next;
5193b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5203b224e63SBarry Smith     }
5213b224e63SBarry Smith   } else {
52297bbdb24SBarry Smith     /* set up the individual PCs */
52397bbdb24SBarry Smith     i    = 0;
5245a9f2f41SSatish Balay     ilink = jac->head;
5255a9f2f41SSatish Balay     while (ilink) {
526519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5273b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
528c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5295a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
53097bbdb24SBarry Smith       i++;
5315a9f2f41SSatish Balay       ilink = ilink->next;
5320971522cSBarry Smith     }
53397bbdb24SBarry Smith 
53497bbdb24SBarry Smith     /* create work vectors for each split */
5351b9fc7fcSBarry Smith     if (!jac->x) {
53697bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5375a9f2f41SSatish Balay       ilink = jac->head;
53897bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
539906ed7ccSBarry Smith         Vec *vl,*vr;
540906ed7ccSBarry Smith 
541906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
542906ed7ccSBarry Smith         ilink->x  = *vr;
543906ed7ccSBarry Smith         ilink->y  = *vl;
544906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
545906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5465a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5475a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5485a9f2f41SSatish Balay         ilink     = ilink->next;
54997bbdb24SBarry Smith       }
5503b224e63SBarry Smith     }
5513b224e63SBarry Smith   }
5523b224e63SBarry Smith 
5533b224e63SBarry Smith 
554d0f46423SBarry Smith   if (!jac->head->sctx) {
5553b224e63SBarry Smith     Vec xtmp;
5563b224e63SBarry Smith 
55779416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5581b9fc7fcSBarry Smith 
5595a9f2f41SSatish Balay     ilink = jac->head;
5601b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5611b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
562704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5635a9f2f41SSatish Balay       ilink = ilink->next;
5641b9fc7fcSBarry Smith     }
5656bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5661b9fc7fcSBarry Smith   }
567c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5680971522cSBarry Smith   PetscFunctionReturn(0);
5690971522cSBarry Smith }
5700971522cSBarry Smith 
5715a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
572ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
573ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5745a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
575ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
576ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
57779416396SBarry Smith 
5780971522cSBarry Smith #undef __FUNCT__
5793b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5803b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5813b224e63SBarry Smith {
5823b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5833b224e63SBarry Smith   PetscErrorCode    ierr;
5843b224e63SBarry Smith   KSP               ksp;
5853b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5863b224e63SBarry Smith 
5873b224e63SBarry Smith   PetscFunctionBegin;
5883b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5893b224e63SBarry Smith 
590c5d2311dSJed Brown   switch (jac->schurfactorization) {
591c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
592a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
593c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
594c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
595c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
596c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
597c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
598c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
599c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
600c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
601c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
602c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
603c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
604c5d2311dSJed Brown       break;
605c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
606a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
607c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
608c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
609c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
610c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
611c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
615c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
616c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
617c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
618c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
619c5d2311dSJed Brown       break;
620c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
621a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
622c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
623c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
624c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
625c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
626c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
627c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
628c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
629c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
630c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
631c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
632c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
633c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
634c5d2311dSJed Brown       break;
635c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
636a04f6461SBarry 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 */
6373b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6383b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6393b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6403b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6413b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6423b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6433b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6443b224e63SBarry Smith 
6453b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6463b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6473b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6483b224e63SBarry Smith 
6493b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6503b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6513b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6523b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6533b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
654c5d2311dSJed Brown   }
6553b224e63SBarry Smith   PetscFunctionReturn(0);
6563b224e63SBarry Smith }
6573b224e63SBarry Smith 
6583b224e63SBarry Smith #undef __FUNCT__
6590971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6600971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6610971522cSBarry Smith {
6620971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6630971522cSBarry Smith   PetscErrorCode    ierr;
6645a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
665af93645bSJed Brown   PetscInt          cnt;
6660971522cSBarry Smith 
6670971522cSBarry Smith   PetscFunctionBegin;
66851f519a2SBarry Smith   CHKMEMQ;
66951f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
67051f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
67151f519a2SBarry Smith 
67279416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6731b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6740971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6755a9f2f41SSatish Balay       while (ilink) {
6765a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6775a9f2f41SSatish Balay         ilink = ilink->next;
6780971522cSBarry Smith       }
6790971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6801b9fc7fcSBarry Smith     } else {
681efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6825a9f2f41SSatish Balay       while (ilink) {
6835a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6845a9f2f41SSatish Balay         ilink = ilink->next;
6851b9fc7fcSBarry Smith       }
6861b9fc7fcSBarry Smith     }
68716913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
68879416396SBarry Smith     if (!jac->w1) {
68979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
69079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
69179416396SBarry Smith     }
692efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6935a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6943e197d65SBarry Smith     cnt = 1;
6955a9f2f41SSatish Balay     while (ilink->next) {
6965a9f2f41SSatish Balay       ilink = ilink->next;
6973e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6983e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6993e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7003e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7013e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7023e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7033e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7043e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7053e197d65SBarry Smith     }
70651f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
70711755939SBarry Smith       cnt -= 2;
70851f519a2SBarry Smith       while (ilink->previous) {
70951f519a2SBarry Smith         ilink = ilink->previous;
71011755939SBarry Smith         /* compute the residual only over the part of the vector needed */
71111755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
71211755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
71311755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71411755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71511755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
71611755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
71711755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
71851f519a2SBarry Smith       }
71911755939SBarry Smith     }
72065e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
72151f519a2SBarry Smith   CHKMEMQ;
7220971522cSBarry Smith   PetscFunctionReturn(0);
7230971522cSBarry Smith }
7240971522cSBarry Smith 
725421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
726ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
727ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
728421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
729ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
730ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
731421e10b8SBarry Smith 
732421e10b8SBarry Smith #undef __FUNCT__
7338c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
734421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
735421e10b8SBarry Smith {
736421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
737421e10b8SBarry Smith   PetscErrorCode    ierr;
738421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
739421e10b8SBarry Smith 
740421e10b8SBarry Smith   PetscFunctionBegin;
741421e10b8SBarry Smith   CHKMEMQ;
742421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
743421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
744421e10b8SBarry Smith 
745421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
746421e10b8SBarry Smith     if (jac->defaultsplit) {
747421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
748421e10b8SBarry Smith       while (ilink) {
749421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
750421e10b8SBarry Smith 	ilink = ilink->next;
751421e10b8SBarry Smith       }
752421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
753421e10b8SBarry Smith     } else {
754421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
755421e10b8SBarry Smith       while (ilink) {
756421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
757421e10b8SBarry Smith 	ilink = ilink->next;
758421e10b8SBarry Smith       }
759421e10b8SBarry Smith     }
760421e10b8SBarry Smith   } else {
761421e10b8SBarry Smith     if (!jac->w1) {
762421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
763421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
764421e10b8SBarry Smith     }
765421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
766421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
767421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
768421e10b8SBarry Smith       while (ilink->next) {
769421e10b8SBarry Smith         ilink = ilink->next;
7709989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
771421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
772421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
773421e10b8SBarry Smith       }
774421e10b8SBarry Smith       while (ilink->previous) {
775421e10b8SBarry Smith         ilink = ilink->previous;
7769989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
777421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
778421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
779421e10b8SBarry Smith       }
780421e10b8SBarry Smith     } else {
781421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
782421e10b8SBarry Smith 	ilink = ilink->next;
783421e10b8SBarry Smith       }
784421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
785421e10b8SBarry Smith       while (ilink->previous) {
786421e10b8SBarry Smith 	ilink = ilink->previous;
7879989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
788421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
789421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
790421e10b8SBarry Smith       }
791421e10b8SBarry Smith     }
792421e10b8SBarry Smith   }
793421e10b8SBarry Smith   CHKMEMQ;
794421e10b8SBarry Smith   PetscFunctionReturn(0);
795421e10b8SBarry Smith }
796421e10b8SBarry Smith 
7970971522cSBarry Smith #undef __FUNCT__
798574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
799574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8000971522cSBarry Smith {
8010971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8020971522cSBarry Smith   PetscErrorCode    ierr;
8035a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8040971522cSBarry Smith 
8050971522cSBarry Smith   PetscFunctionBegin;
8065a9f2f41SSatish Balay   while (ilink) {
807574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
808fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
809fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
810fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
811fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
8125a9f2f41SSatish Balay     next = ilink->next;
8135a9f2f41SSatish Balay     ilink = next;
8140971522cSBarry Smith   }
81505b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
816574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
817574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
818574deadeSBarry Smith   } else if (jac->mat) {
819574deadeSBarry Smith     jac->mat = PETSC_NULL;
820574deadeSBarry Smith   }
82197bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
82268dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8236bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8246bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8256bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8266bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8276bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8286bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8296bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
83063ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
831574deadeSBarry Smith   PetscFunctionReturn(0);
832574deadeSBarry Smith }
833574deadeSBarry Smith 
834574deadeSBarry Smith #undef __FUNCT__
835574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
836574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
837574deadeSBarry Smith {
838574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
839574deadeSBarry Smith   PetscErrorCode    ierr;
840574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
841574deadeSBarry Smith 
842574deadeSBarry Smith   PetscFunctionBegin;
843574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
844574deadeSBarry Smith   while (ilink) {
8456bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
846574deadeSBarry Smith     next = ilink->next;
847574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
848574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
849574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
850574deadeSBarry Smith     ilink = next;
851574deadeSBarry Smith   }
852574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
853c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
854d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
855d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
856d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
857d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
858d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
8590971522cSBarry Smith   PetscFunctionReturn(0);
8600971522cSBarry Smith }
8610971522cSBarry Smith 
8620971522cSBarry Smith #undef __FUNCT__
8630971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8640971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8650971522cSBarry Smith {
8661b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8676c924f48SJed Brown   PetscInt        bs;
868bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8699dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8703b224e63SBarry Smith   PCCompositeType ctype;
8711b9fc7fcSBarry Smith 
8720971522cSBarry Smith   PetscFunctionBegin;
8731b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
874acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
87551f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
87651f519a2SBarry Smith   if (flg) {
87751f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
87851f519a2SBarry Smith   }
879704ba839SBarry Smith 
880435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
881c0adfefeSBarry Smith   if (stokes) {
882c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
883c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
884c0adfefeSBarry Smith   }
885c0adfefeSBarry Smith 
8863b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8873b224e63SBarry Smith   if (flg) {
8883b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8893b224e63SBarry Smith   }
890d32f9abdSBarry Smith 
891c30613efSMatthew Knepley   /* Only setup fields once */
892c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
893d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
894d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8956c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8966c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
897d32f9abdSBarry Smith   }
898c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
899c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
900084e4875SJed 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);
901c5d2311dSJed Brown   }
9021b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9030971522cSBarry Smith   PetscFunctionReturn(0);
9040971522cSBarry Smith }
9050971522cSBarry Smith 
9060971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9070971522cSBarry Smith 
9080971522cSBarry Smith EXTERN_C_BEGIN
9090971522cSBarry Smith #undef __FUNCT__
9100971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9117087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
9120971522cSBarry Smith {
91397bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9140971522cSBarry Smith   PetscErrorCode    ierr;
9155a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
91669a612a9SBarry Smith   char              prefix[128];
91751f519a2SBarry Smith   PetscInt          i;
9180971522cSBarry Smith 
9190971522cSBarry Smith   PetscFunctionBegin;
9206c924f48SJed Brown   if (jac->splitdefined) {
9216c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9226c924f48SJed Brown     PetscFunctionReturn(0);
9236c924f48SJed Brown   }
92451f519a2SBarry Smith   for (i=0; i<n; i++) {
925e32f2f54SBarry 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);
926e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
92751f519a2SBarry Smith   }
928704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
929a04f6461SBarry Smith   if (splitname) {
930db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
931a04f6461SBarry Smith   } else {
932a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
933a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
934a04f6461SBarry Smith   }
935704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9365a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9375a9f2f41SSatish Balay   ilink->nfields = n;
9385a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9397adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9401cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9415a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9429005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
94369a612a9SBarry Smith 
944a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9455a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9460971522cSBarry Smith 
9470971522cSBarry Smith   if (!next) {
9485a9f2f41SSatish Balay     jac->head       = ilink;
94951f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9500971522cSBarry Smith   } else {
9510971522cSBarry Smith     while (next->next) {
9520971522cSBarry Smith       next = next->next;
9530971522cSBarry Smith     }
9545a9f2f41SSatish Balay     next->next      = ilink;
95551f519a2SBarry Smith     ilink->previous = next;
9560971522cSBarry Smith   }
9570971522cSBarry Smith   jac->nsplits++;
9580971522cSBarry Smith   PetscFunctionReturn(0);
9590971522cSBarry Smith }
9600971522cSBarry Smith EXTERN_C_END
9610971522cSBarry Smith 
962e69d4d44SBarry Smith EXTERN_C_BEGIN
963e69d4d44SBarry Smith #undef __FUNCT__
964e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9657087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
966e69d4d44SBarry Smith {
967e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
968e69d4d44SBarry Smith   PetscErrorCode ierr;
969e69d4d44SBarry Smith 
970e69d4d44SBarry Smith   PetscFunctionBegin;
971519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
972e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
973e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
97413e0d083SBarry Smith   if (n) *n = jac->nsplits;
975e69d4d44SBarry Smith   PetscFunctionReturn(0);
976e69d4d44SBarry Smith }
977e69d4d44SBarry Smith EXTERN_C_END
9780971522cSBarry Smith 
9790971522cSBarry Smith EXTERN_C_BEGIN
9800971522cSBarry Smith #undef __FUNCT__
98169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9827087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9830971522cSBarry Smith {
9840971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9850971522cSBarry Smith   PetscErrorCode    ierr;
9860971522cSBarry Smith   PetscInt          cnt = 0;
9875a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9880971522cSBarry Smith 
9890971522cSBarry Smith   PetscFunctionBegin;
9905d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9915a9f2f41SSatish Balay   while (ilink) {
9925a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9935a9f2f41SSatish Balay     ilink = ilink->next;
9940971522cSBarry Smith   }
9955d480477SMatthew 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);
99613e0d083SBarry Smith   if (n) *n = jac->nsplits;
9970971522cSBarry Smith   PetscFunctionReturn(0);
9980971522cSBarry Smith }
9990971522cSBarry Smith EXTERN_C_END
10000971522cSBarry Smith 
1001704ba839SBarry Smith EXTERN_C_BEGIN
1002704ba839SBarry Smith #undef __FUNCT__
1003704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10047087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1005704ba839SBarry Smith {
1006704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1007704ba839SBarry Smith   PetscErrorCode    ierr;
1008704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1009704ba839SBarry Smith   char              prefix[128];
1010704ba839SBarry Smith 
1011704ba839SBarry Smith   PetscFunctionBegin;
10126c924f48SJed Brown   if (jac->splitdefined) {
10136c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10146c924f48SJed Brown     PetscFunctionReturn(0);
10156c924f48SJed Brown   }
101616913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1017a04f6461SBarry Smith   if (splitname) {
1018db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1019a04f6461SBarry Smith   } else {
1020a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1021d88c663fSJed Brown     ierr = PetscSNPrintf(ilink->splitname,3,"%D",jac->nsplits);CHKERRQ(ierr);
1022a04f6461SBarry Smith   }
10231ab39975SBarry Smith   ilink->is      = is;
10241ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1025704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1026704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10271cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1028704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10299005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1030704ba839SBarry Smith 
1031a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1032704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1033704ba839SBarry Smith 
1034704ba839SBarry Smith   if (!next) {
1035704ba839SBarry Smith     jac->head       = ilink;
1036704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1037704ba839SBarry Smith   } else {
1038704ba839SBarry Smith     while (next->next) {
1039704ba839SBarry Smith       next = next->next;
1040704ba839SBarry Smith     }
1041704ba839SBarry Smith     next->next      = ilink;
1042704ba839SBarry Smith     ilink->previous = next;
1043704ba839SBarry Smith   }
1044704ba839SBarry Smith   jac->nsplits++;
1045704ba839SBarry Smith 
1046704ba839SBarry Smith   PetscFunctionReturn(0);
1047704ba839SBarry Smith }
1048704ba839SBarry Smith EXTERN_C_END
1049704ba839SBarry Smith 
10500971522cSBarry Smith #undef __FUNCT__
10510971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10520971522cSBarry Smith /*@
10530971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10540971522cSBarry Smith 
1055ad4df100SBarry Smith     Logically Collective on PC
10560971522cSBarry Smith 
10570971522cSBarry Smith     Input Parameters:
10580971522cSBarry Smith +   pc  - the preconditioner context
1059a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10600971522cSBarry Smith .   n - the number of fields in this split
1061db4c96c1SJed Brown -   fields - the fields in this split
10620971522cSBarry Smith 
10630971522cSBarry Smith     Level: intermediate
10640971522cSBarry Smith 
1065d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1066d32f9abdSBarry Smith 
1067d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1068d32f9abdSBarry 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
1069d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1070d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1071d32f9abdSBarry Smith 
1072db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1073db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1074db4c96c1SJed Brown 
1075d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
10760971522cSBarry Smith 
10770971522cSBarry Smith @*/
10787087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10790971522cSBarry Smith {
10804ac538c5SBarry Smith   PetscErrorCode ierr;
10810971522cSBarry Smith 
10820971522cSBarry Smith   PetscFunctionBegin;
10830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1084db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1085db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1086db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10874ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10880971522cSBarry Smith   PetscFunctionReturn(0);
10890971522cSBarry Smith }
10900971522cSBarry Smith 
10910971522cSBarry Smith #undef __FUNCT__
1092704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1093704ba839SBarry Smith /*@
1094704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1095704ba839SBarry Smith 
1096ad4df100SBarry Smith     Logically Collective on PC
1097704ba839SBarry Smith 
1098704ba839SBarry Smith     Input Parameters:
1099704ba839SBarry Smith +   pc  - the preconditioner context
1100a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1101db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1102704ba839SBarry Smith 
1103d32f9abdSBarry Smith 
1104a6ffb8dbSJed Brown     Notes:
1105a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1106a6ffb8dbSJed Brown 
1107db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1108db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1109d32f9abdSBarry Smith 
1110704ba839SBarry Smith     Level: intermediate
1111704ba839SBarry Smith 
1112704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1113704ba839SBarry Smith 
1114704ba839SBarry Smith @*/
11157087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1116704ba839SBarry Smith {
11174ac538c5SBarry Smith   PetscErrorCode ierr;
1118704ba839SBarry Smith 
1119704ba839SBarry Smith   PetscFunctionBegin;
11200700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1121d88c663fSJed Brown   if (splitname) PetscValidCharPointer(splitname,2);
1122db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11234ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1124704ba839SBarry Smith   PetscFunctionReturn(0);
1125704ba839SBarry Smith }
1126704ba839SBarry Smith 
1127704ba839SBarry Smith #undef __FUNCT__
112857a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
112957a9adfeSMatthew G Knepley /*@
113057a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
113157a9adfeSMatthew G Knepley 
113257a9adfeSMatthew G Knepley     Logically Collective on PC
113357a9adfeSMatthew G Knepley 
113457a9adfeSMatthew G Knepley     Input Parameters:
113557a9adfeSMatthew G Knepley +   pc  - the preconditioner context
113657a9adfeSMatthew G Knepley -   splitname - name of this split
113757a9adfeSMatthew G Knepley 
113857a9adfeSMatthew G Knepley     Output Parameter:
113957a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
114057a9adfeSMatthew G Knepley 
114157a9adfeSMatthew G Knepley     Level: intermediate
114257a9adfeSMatthew G Knepley 
114357a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
114457a9adfeSMatthew G Knepley 
114557a9adfeSMatthew G Knepley @*/
114657a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
114757a9adfeSMatthew G Knepley {
114857a9adfeSMatthew G Knepley   PetscErrorCode ierr;
114957a9adfeSMatthew G Knepley 
115057a9adfeSMatthew G Knepley   PetscFunctionBegin;
115157a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115257a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
115357a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
115457a9adfeSMatthew G Knepley   {
115557a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
115657a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
115757a9adfeSMatthew G Knepley     PetscBool         found;
115857a9adfeSMatthew G Knepley 
115957a9adfeSMatthew G Knepley     *is = PETSC_NULL;
116057a9adfeSMatthew G Knepley     while(ilink) {
116157a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
116257a9adfeSMatthew G Knepley       if (found) {
116357a9adfeSMatthew G Knepley         *is = ilink->is;
116457a9adfeSMatthew G Knepley         break;
116557a9adfeSMatthew G Knepley       }
116657a9adfeSMatthew G Knepley       ilink = ilink->next;
116757a9adfeSMatthew G Knepley     }
116857a9adfeSMatthew G Knepley   }
116957a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
117057a9adfeSMatthew G Knepley }
117157a9adfeSMatthew G Knepley 
117257a9adfeSMatthew G Knepley #undef __FUNCT__
117351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
117451f519a2SBarry Smith /*@
117551f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
117651f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
117751f519a2SBarry Smith 
1178ad4df100SBarry Smith     Logically Collective on PC
117951f519a2SBarry Smith 
118051f519a2SBarry Smith     Input Parameters:
118151f519a2SBarry Smith +   pc  - the preconditioner context
118251f519a2SBarry Smith -   bs - the block size
118351f519a2SBarry Smith 
118451f519a2SBarry Smith     Level: intermediate
118551f519a2SBarry Smith 
118651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
118751f519a2SBarry Smith 
118851f519a2SBarry Smith @*/
11897087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
119051f519a2SBarry Smith {
11914ac538c5SBarry Smith   PetscErrorCode ierr;
119251f519a2SBarry Smith 
119351f519a2SBarry Smith   PetscFunctionBegin;
11940700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1195c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
11964ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
119751f519a2SBarry Smith   PetscFunctionReturn(0);
119851f519a2SBarry Smith }
119951f519a2SBarry Smith 
120051f519a2SBarry Smith #undef __FUNCT__
120169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12020971522cSBarry Smith /*@C
120369a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12040971522cSBarry Smith 
120569a612a9SBarry Smith    Collective on KSP
12060971522cSBarry Smith 
12070971522cSBarry Smith    Input Parameter:
12080971522cSBarry Smith .  pc - the preconditioner context
12090971522cSBarry Smith 
12100971522cSBarry Smith    Output Parameters:
121113e0d083SBarry Smith +  n - the number of splits
121269a612a9SBarry Smith -  pc - the array of KSP contexts
12130971522cSBarry Smith 
12140971522cSBarry Smith    Note:
1215d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1216d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12170971522cSBarry Smith 
121869a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12190971522cSBarry Smith 
12200971522cSBarry Smith    Level: advanced
12210971522cSBarry Smith 
12220971522cSBarry Smith .seealso: PCFIELDSPLIT
12230971522cSBarry Smith @*/
12247087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12250971522cSBarry Smith {
12264ac538c5SBarry Smith   PetscErrorCode ierr;
12270971522cSBarry Smith 
12280971522cSBarry Smith   PetscFunctionBegin;
12290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
123013e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12314ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12320971522cSBarry Smith   PetscFunctionReturn(0);
12330971522cSBarry Smith }
12340971522cSBarry Smith 
1235e69d4d44SBarry Smith #undef __FUNCT__
1236e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1237e69d4d44SBarry Smith /*@
1238e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1239a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1240e69d4d44SBarry Smith 
1241e69d4d44SBarry Smith     Collective on PC
1242e69d4d44SBarry Smith 
1243e69d4d44SBarry Smith     Input Parameters:
1244e69d4d44SBarry Smith +   pc  - the preconditioner context
1245fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1246084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1247084e4875SJed Brown 
1248e69d4d44SBarry Smith     Options Database:
1249084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1250e69d4d44SBarry Smith 
1251fd1303e9SJungho Lee     Notes:
1252fd1303e9SJungho Lee $    If ptype is
1253fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1254fd1303e9SJungho Lee $             to this function).
1255fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1256fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1257fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1258fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1259fd1303e9SJungho Lee $             preconditioner
1260fd1303e9SJungho Lee 
1261fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1262fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1263fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1264fd1303e9SJungho Lee 
1265fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1266fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1267fd1303e9SJungho Lee 
1268e69d4d44SBarry Smith     Level: intermediate
1269e69d4d44SBarry Smith 
1270fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1271e69d4d44SBarry Smith 
1272e69d4d44SBarry Smith @*/
12737087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1274e69d4d44SBarry Smith {
12754ac538c5SBarry Smith   PetscErrorCode ierr;
1276e69d4d44SBarry Smith 
1277e69d4d44SBarry Smith   PetscFunctionBegin;
12780700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12794ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1280e69d4d44SBarry Smith   PetscFunctionReturn(0);
1281e69d4d44SBarry Smith }
1282e69d4d44SBarry Smith 
1283e69d4d44SBarry Smith EXTERN_C_BEGIN
1284e69d4d44SBarry Smith #undef __FUNCT__
1285e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
12867087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1287e69d4d44SBarry Smith {
1288e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1289084e4875SJed Brown   PetscErrorCode  ierr;
1290e69d4d44SBarry Smith 
1291e69d4d44SBarry Smith   PetscFunctionBegin;
1292084e4875SJed Brown   jac->schurpre = ptype;
1293084e4875SJed Brown   if (pre) {
12946bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1295084e4875SJed Brown     jac->schur_user = pre;
1296084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1297084e4875SJed Brown   }
1298e69d4d44SBarry Smith   PetscFunctionReturn(0);
1299e69d4d44SBarry Smith }
1300e69d4d44SBarry Smith EXTERN_C_END
1301e69d4d44SBarry Smith 
130230ad9308SMatthew Knepley #undef __FUNCT__
130330ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
130430ad9308SMatthew Knepley /*@C
13058c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
130630ad9308SMatthew Knepley 
130730ad9308SMatthew Knepley    Collective on KSP
130830ad9308SMatthew Knepley 
130930ad9308SMatthew Knepley    Input Parameter:
131030ad9308SMatthew Knepley .  pc - the preconditioner context
131130ad9308SMatthew Knepley 
131230ad9308SMatthew Knepley    Output Parameters:
1313a04f6461SBarry Smith +  A00 - the (0,0) block
1314a04f6461SBarry Smith .  A01 - the (0,1) block
1315a04f6461SBarry Smith .  A10 - the (1,0) block
1316a04f6461SBarry Smith -  A11 - the (1,1) block
131730ad9308SMatthew Knepley 
131830ad9308SMatthew Knepley    Level: advanced
131930ad9308SMatthew Knepley 
132030ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
132130ad9308SMatthew Knepley @*/
1322a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
132330ad9308SMatthew Knepley {
132430ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
132530ad9308SMatthew Knepley 
132630ad9308SMatthew Knepley   PetscFunctionBegin;
13270700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1328c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1329a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1330a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1331a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1332a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
133330ad9308SMatthew Knepley   PetscFunctionReturn(0);
133430ad9308SMatthew Knepley }
133530ad9308SMatthew Knepley 
133679416396SBarry Smith EXTERN_C_BEGIN
133779416396SBarry Smith #undef __FUNCT__
133879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
13397087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
134079416396SBarry Smith {
134179416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1342e69d4d44SBarry Smith   PetscErrorCode ierr;
134379416396SBarry Smith 
134479416396SBarry Smith   PetscFunctionBegin;
134579416396SBarry Smith   jac->type = type;
13463b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
13473b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
13483b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1349e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1350e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1351e69d4d44SBarry Smith 
13523b224e63SBarry Smith   } else {
13533b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
13543b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1355e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13569e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13573b224e63SBarry Smith   }
135879416396SBarry Smith   PetscFunctionReturn(0);
135979416396SBarry Smith }
136079416396SBarry Smith EXTERN_C_END
136179416396SBarry Smith 
136251f519a2SBarry Smith EXTERN_C_BEGIN
136351f519a2SBarry Smith #undef __FUNCT__
136451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13657087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
136651f519a2SBarry Smith {
136751f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
136851f519a2SBarry Smith 
136951f519a2SBarry Smith   PetscFunctionBegin;
137065e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
137165e19b50SBarry 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);
137251f519a2SBarry Smith   jac->bs = bs;
137351f519a2SBarry Smith   PetscFunctionReturn(0);
137451f519a2SBarry Smith }
137551f519a2SBarry Smith EXTERN_C_END
137651f519a2SBarry Smith 
137779416396SBarry Smith #undef __FUNCT__
137879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1379bc08b0f1SBarry Smith /*@
138079416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
138179416396SBarry Smith 
138279416396SBarry Smith    Collective on PC
138379416396SBarry Smith 
138479416396SBarry Smith    Input Parameter:
138579416396SBarry Smith .  pc - the preconditioner context
138681540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
138779416396SBarry Smith 
138879416396SBarry Smith    Options Database Key:
1389a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
139079416396SBarry Smith 
139179416396SBarry Smith    Level: Developer
139279416396SBarry Smith 
139379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
139479416396SBarry Smith 
139579416396SBarry Smith .seealso: PCCompositeSetType()
139679416396SBarry Smith 
139779416396SBarry Smith @*/
13987087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
139979416396SBarry Smith {
14004ac538c5SBarry Smith   PetscErrorCode ierr;
140179416396SBarry Smith 
140279416396SBarry Smith   PetscFunctionBegin;
14030700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14044ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
140579416396SBarry Smith   PetscFunctionReturn(0);
140679416396SBarry Smith }
140779416396SBarry Smith 
14080971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
14090971522cSBarry Smith /*MC
1410a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1411a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
14120971522cSBarry Smith 
1413edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1414edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
14150971522cSBarry Smith 
1416a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
141769a612a9SBarry Smith          and set the options directly on the resulting KSP object
14180971522cSBarry Smith 
14190971522cSBarry Smith    Level: intermediate
14200971522cSBarry Smith 
142179416396SBarry Smith    Options Database Keys:
142281540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
142381540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
142481540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
142581540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
142681540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1427e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1428435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
142943890ad2SJed Brown .   -fieldsplit_NAME_ksp_* - control inner linear solver, NAME is a sequential integer if unspecified, otherwise use name provided in PCFieldSplitSetIS() or the name associated with the field in the DM passed to PCSetDM() (or a higher level function)
143043890ad2SJed Brown -   -fieldsplit_NAME_pc_* - control inner preconditioner, NAME has same semantics as above
143179416396SBarry Smith 
143243890ad2SJed Brown    Notes:
143343890ad2SJed Brown     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1434d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1435d32f9abdSBarry Smith 
1436d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1437d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1438d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1439d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1440d32f9abdSBarry Smith 
1441a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1442a04f6461SBarry Smith                                                     ( A10 A11 )
1443e69d4d44SBarry Smith      the preconditioner is
1444a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1445a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1446a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1447a04f6461SBarry 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).
1448a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1449edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1450a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
14517e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1452e69d4d44SBarry Smith 
1453edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1454edf189efSBarry Smith      is used automatically for a second block.
1455edf189efSBarry Smith 
1456ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1457ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1458ff218e97SBarry Smith 
1459ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1460ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1461ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
14620716a85fSBarry Smith 
1463a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
14640971522cSBarry Smith 
14657e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1466e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
14670971522cSBarry Smith M*/
14680971522cSBarry Smith 
14690971522cSBarry Smith EXTERN_C_BEGIN
14700971522cSBarry Smith #undef __FUNCT__
14710971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
14727087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
14730971522cSBarry Smith {
14740971522cSBarry Smith   PetscErrorCode ierr;
14750971522cSBarry Smith   PC_FieldSplit  *jac;
14760971522cSBarry Smith 
14770971522cSBarry Smith   PetscFunctionBegin;
147838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
14790971522cSBarry Smith   jac->bs        = -1;
14800971522cSBarry Smith   jac->nsplits   = 0;
14813e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1482e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1483c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
148451f519a2SBarry Smith 
14850971522cSBarry Smith   pc->data     = (void*)jac;
14860971522cSBarry Smith 
14870971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1488421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
14890971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1490574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
14910971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
14920971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
14930971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
14940971522cSBarry Smith   pc->ops->applyrichardson   = 0;
14950971522cSBarry Smith 
149669a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
149769a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14980971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
14990971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1500704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1501704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
150279416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
150379416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
150451f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
150551f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
15060971522cSBarry Smith   PetscFunctionReturn(0);
15070971522cSBarry Smith }
15080971522cSBarry Smith EXTERN_C_END
15090971522cSBarry Smith 
15100971522cSBarry Smith 
1511a541d17aSBarry Smith 
1512