xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 5668aaf46f6adb7d49c9f2a9153c8d067db11cfe)
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;
215d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
221b9fc7fcSBarry Smith   VecScatter        sctx;
235d4c12cdSJungho Lee   IS                is,is_col;
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;
1785d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
1795d4c12cdSJungho Lee   PetscBool      flg,flg_col;
1805d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
1816c924f48SJed Brown 
1826c924f48SJed Brown   PetscFunctionBegin;
1836c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1845d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
1856c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1866c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1876c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1885d4c12cdSJungho Lee     ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
1896c924f48SJed Brown     nfields = jac->bs;
19029499fbbSJungho Lee     nfields_col = jac->bs;
1916c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1925d4c12cdSJungho Lee     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
1936c924f48SJed Brown     if (!flg) break;
1945d4c12cdSJungho Lee     else if (flg && !flg_col) {
1955d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1965d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
1975d4c12cdSJungho Lee     }
1985d4c12cdSJungho Lee     else {
1995d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2005d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
2015d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
2025d4c12cdSJungho Lee     }
2036c924f48SJed Brown   }
2046c924f48SJed Brown   if (i > 0) {
2056c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
2066c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2076c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2086c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2096c924f48SJed Brown   }
2106c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2115d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2126c924f48SJed Brown   PetscFunctionReturn(0);
2136c924f48SJed Brown }
2146c924f48SJed Brown 
2156c924f48SJed Brown #undef __FUNCT__
21669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
21769a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2180971522cSBarry Smith {
2190971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2200971522cSBarry Smith   PetscErrorCode    ierr;
2215a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2226ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2236c924f48SJed Brown   PetscInt          i;
2240971522cSBarry Smith 
2250971522cSBarry Smith   PetscFunctionBegin;
226d32f9abdSBarry Smith   if (!ilink) {
227d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
228d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2297b62db95SJungho Lee       PetscInt     numFields, f;
2300784a22cSJed Brown       char         **fieldNames;
2317b62db95SJungho Lee       IS          *fields;
2324d343eeaSMatthew G Knepley       PetscBool    dmcomposite;
2337b62db95SJungho Lee 
2347b62db95SJungho Lee       ierr = DMCreateFieldIS(pc->dm, &numFields, &fieldNames, &fields);CHKERRQ(ierr);
2357b62db95SJungho Lee       for(f = 0; f < numFields; ++f) {
2367b62db95SJungho Lee         ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
2377b62db95SJungho Lee         ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
2387b62db95SJungho Lee         ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
2397b62db95SJungho Lee       }
2407b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
2417b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
2427b62db95SJungho Lee 
2434d343eeaSMatthew G Knepley       ierr = PetscTypeCompare((PetscObject) pc->dm, DMCOMPOSITE, &dmcomposite);CHKERRQ(ierr);
2444d343eeaSMatthew G Knepley       if (dmcomposite) {
2457b62db95SJungho Lee         DM      *dms;
2464d343eeaSMatthew G Knepley         PetscInt nDM;
2477b62db95SJungho Lee 
2488b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2498b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm, &nDM);CHKERRQ(ierr);
2507b62db95SJungho Lee         ierr = PetscMalloc(nDM*sizeof(DM), &dms);CHKERRQ(ierr);
2517b62db95SJungho Lee         ierr = DMCompositeGetEntriesArray(pc->dm, dms);CHKERRQ(ierr);
2527b62db95SJungho Lee         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2537b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
2547b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
2552fa5ba8aSJed Brown         }
2567b62db95SJungho Lee         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);
2915d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
29279416396SBarry Smith           }
2935d4c12cdSJungho Lee           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);
3067b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,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;
3187b62db95SJungho Lee         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);
3227b62db95SJungho Lee         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3237b62db95SJungho Lee         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
324d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3257b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3267b62db95SJungho Lee           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);
3317b62db95SJungho Lee         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 
3507b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
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;
3635f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
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 
3765d4c12cdSJungho Lee     /* 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);
3845f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
385704ba839SBarry Smith       } else if (!ilink->is) {
386ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3875f4ab4e1SJungho Lee           PetscInt   *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
3885a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3895f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
3901b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3911b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3921b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
3935f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
39497bbdb24SBarry Smith             }
39597bbdb24SBarry Smith           }
396d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
3975f4ab4e1SJungho Lee           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
398ccb205f8SBarry Smith         } else {
399704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
4005f4ab4e1SJungho Lee           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
401ccb205f8SBarry Smith        }
4023e197d65SBarry Smith       }
403704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
4045f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
4055f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
406704ba839SBarry Smith       ilink = ilink->next;
4071b9fc7fcSBarry Smith     }
4081b9fc7fcSBarry Smith   }
4091b9fc7fcSBarry Smith 
410704ba839SBarry Smith   ilink  = jac->head;
41197bbdb24SBarry Smith   if (!jac->pmat) {
412cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
413cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4145f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
415704ba839SBarry Smith       ilink = ilink->next;
416cf502942SBarry Smith     }
41797bbdb24SBarry Smith   } else {
418cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4195f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
420704ba839SBarry Smith       ilink = ilink->next;
421cf502942SBarry Smith     }
42297bbdb24SBarry Smith   }
423519d70e2SJed Brown   if (jac->realdiagonal) {
424519d70e2SJed Brown     ilink = jac->head;
425519d70e2SJed Brown     if (!jac->mat) {
426519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
427519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4285f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
429519d70e2SJed Brown         ilink = ilink->next;
430519d70e2SJed Brown       }
431519d70e2SJed Brown     } else {
432519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4335f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
434519d70e2SJed Brown         ilink = ilink->next;
435519d70e2SJed Brown       }
436519d70e2SJed Brown     }
437519d70e2SJed Brown   } else {
438519d70e2SJed Brown     jac->mat = jac->pmat;
439519d70e2SJed Brown   }
44097bbdb24SBarry Smith 
4416c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
44268dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
44368dd23aaSBarry Smith     ilink  = jac->head;
44468dd23aaSBarry Smith     if (!jac->Afield) {
44568dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
44668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4474aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
44868dd23aaSBarry Smith         ilink = ilink->next;
44968dd23aaSBarry Smith       }
45068dd23aaSBarry Smith     } else {
45168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4524aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
45368dd23aaSBarry Smith         ilink = ilink->next;
45468dd23aaSBarry Smith       }
45568dd23aaSBarry Smith     }
45668dd23aaSBarry Smith   }
45768dd23aaSBarry Smith 
4583b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4593b224e63SBarry Smith     IS       ccis;
4604aa3045dSJed Brown     PetscInt rstart,rend;
461093c86ffSJed Brown     char     lscname[256];
462093c86ffSJed Brown     PetscObject LSC_L;
463e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
46468dd23aaSBarry Smith 
465e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
466e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
467e6cab6aaSJed Brown 
4683b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4693b224e63SBarry Smith     if (jac->schur) {
4703b224e63SBarry Smith       ilink = jac->head;
47149bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4724aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
473fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4743b224e63SBarry Smith       ilink = ilink->next;
47549bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4764aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
477fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
478519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
479084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4803b224e63SBarry Smith 
4813b224e63SBarry Smith      } else {
4821cee3971SBarry Smith       KSP ksp;
4836c924f48SJed Brown       char schurprefix[256];
4843b224e63SBarry Smith 
485a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4863b224e63SBarry Smith       ilink = jac->head;
48749bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4884aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
489fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4903b224e63SBarry Smith       ilink = ilink->next;
49149bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4924aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
493fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4947d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
495519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
496a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4971cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
498e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
499a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
500a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
50120b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
50220b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
5037b62db95SJungho Lee       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
5043b224e63SBarry Smith 
5053b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
5069005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
5071cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
508084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
509084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
510084e4875SJed Brown         PC pc;
511cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
512084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
513084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
514e69d4d44SBarry Smith       }
5156c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5166c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5173b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
51820b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
51920b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5203b224e63SBarry Smith 
5213b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5223b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5233b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5243b224e63SBarry Smith       ilink = jac->head;
5253b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5263b224e63SBarry Smith       ilink = ilink->next;
5273b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5283b224e63SBarry Smith     }
529093c86ffSJed Brown 
530093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
531093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
532093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
533093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
534093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
535093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
536093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
537093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
538093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
5393b224e63SBarry Smith   } else {
54097bbdb24SBarry Smith     /* set up the individual PCs */
54197bbdb24SBarry Smith     i    = 0;
5425a9f2f41SSatish Balay     ilink = jac->head;
5435a9f2f41SSatish Balay     while (ilink) {
544519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5453b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
546c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5475a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
54897bbdb24SBarry Smith       i++;
5495a9f2f41SSatish Balay       ilink = ilink->next;
5500971522cSBarry Smith     }
55197bbdb24SBarry Smith 
55297bbdb24SBarry Smith     /* create work vectors for each split */
5531b9fc7fcSBarry Smith     if (!jac->x) {
55497bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5555a9f2f41SSatish Balay       ilink = jac->head;
55697bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
557906ed7ccSBarry Smith         Vec *vl,*vr;
558906ed7ccSBarry Smith 
559906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
560906ed7ccSBarry Smith         ilink->x  = *vr;
561906ed7ccSBarry Smith         ilink->y  = *vl;
562906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
563906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5645a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5655a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5665a9f2f41SSatish Balay         ilink     = ilink->next;
56797bbdb24SBarry Smith       }
5683b224e63SBarry Smith     }
5693b224e63SBarry Smith   }
5703b224e63SBarry Smith 
5713b224e63SBarry Smith 
572d0f46423SBarry Smith   if (!jac->head->sctx) {
5733b224e63SBarry Smith     Vec xtmp;
5743b224e63SBarry Smith 
57579416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5761b9fc7fcSBarry Smith 
5775a9f2f41SSatish Balay     ilink = jac->head;
5781b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5791b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
580704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5815a9f2f41SSatish Balay       ilink = ilink->next;
5821b9fc7fcSBarry Smith     }
5836bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5841b9fc7fcSBarry Smith   }
585c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5860971522cSBarry Smith   PetscFunctionReturn(0);
5870971522cSBarry Smith }
5880971522cSBarry Smith 
5895a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
590ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
591ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5925a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
593ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
594ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
59579416396SBarry Smith 
5960971522cSBarry Smith #undef __FUNCT__
5973b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5983b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5993b224e63SBarry Smith {
6003b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6013b224e63SBarry Smith   PetscErrorCode    ierr;
6023b224e63SBarry Smith   KSP               ksp;
6033b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
6043b224e63SBarry Smith 
6053b224e63SBarry Smith   PetscFunctionBegin;
6063b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
6073b224e63SBarry Smith 
608c5d2311dSJed Brown   switch (jac->schurfactorization) {
609c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
610a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
611c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
615c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
616c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
617c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
618c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
619c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
620c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
621c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
622c5d2311dSJed Brown       break;
623c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
624a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
625c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
628c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
629c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
630c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
631c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
632c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
633c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
634c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
635c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
636c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
637c5d2311dSJed Brown       break;
638c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
639a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
640c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
641c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
642c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
643c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
644c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
645c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
646c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
647c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
648c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
649c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
650c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
651c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
652c5d2311dSJed Brown       break;
653c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
654a04f6461SBarry 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 */
6553b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6563b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6573b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6583b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6593b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6603b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6613b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6623b224e63SBarry Smith 
6633b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6643b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6653b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6663b224e63SBarry Smith 
6673b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6683b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6693b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6703b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6713b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
672c5d2311dSJed Brown   }
6733b224e63SBarry Smith   PetscFunctionReturn(0);
6743b224e63SBarry Smith }
6753b224e63SBarry Smith 
6763b224e63SBarry Smith #undef __FUNCT__
6770971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6780971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6790971522cSBarry Smith {
6800971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6810971522cSBarry Smith   PetscErrorCode    ierr;
6825a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
683af93645bSJed Brown   PetscInt          cnt;
6840971522cSBarry Smith 
6850971522cSBarry Smith   PetscFunctionBegin;
68651f519a2SBarry Smith   CHKMEMQ;
68751f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
68851f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
68951f519a2SBarry Smith 
69079416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6911b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6920971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6935a9f2f41SSatish Balay       while (ilink) {
6945a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6955a9f2f41SSatish Balay         ilink = ilink->next;
6960971522cSBarry Smith       }
6970971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6981b9fc7fcSBarry Smith     } else {
699efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
7005a9f2f41SSatish Balay       while (ilink) {
7015a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7025a9f2f41SSatish Balay         ilink = ilink->next;
7031b9fc7fcSBarry Smith       }
7041b9fc7fcSBarry Smith     }
70516913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
70679416396SBarry Smith     if (!jac->w1) {
70779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
70879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
70979416396SBarry Smith     }
710efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
7115a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7123e197d65SBarry Smith     cnt = 1;
7135a9f2f41SSatish Balay     while (ilink->next) {
7145a9f2f41SSatish Balay       ilink = ilink->next;
7153e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7163e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7173e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7183e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7193e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7203e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7213e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7223e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7233e197d65SBarry Smith     }
72451f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
72511755939SBarry Smith       cnt -= 2;
72651f519a2SBarry Smith       while (ilink->previous) {
72751f519a2SBarry Smith         ilink = ilink->previous;
72811755939SBarry Smith         /* compute the residual only over the part of the vector needed */
72911755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
73011755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
73111755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
73211755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
73311755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
73411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
73511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
73651f519a2SBarry Smith       }
73711755939SBarry Smith     }
73865e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
73951f519a2SBarry Smith   CHKMEMQ;
7400971522cSBarry Smith   PetscFunctionReturn(0);
7410971522cSBarry Smith }
7420971522cSBarry Smith 
743421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
744ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
745ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
746421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
747ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
748ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
749421e10b8SBarry Smith 
750421e10b8SBarry Smith #undef __FUNCT__
7518c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
752421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
753421e10b8SBarry Smith {
754421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
755421e10b8SBarry Smith   PetscErrorCode    ierr;
756421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
757421e10b8SBarry Smith 
758421e10b8SBarry Smith   PetscFunctionBegin;
759421e10b8SBarry Smith   CHKMEMQ;
760421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
761421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
762421e10b8SBarry Smith 
763421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
764421e10b8SBarry Smith     if (jac->defaultsplit) {
765421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
766421e10b8SBarry Smith       while (ilink) {
767421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
768421e10b8SBarry Smith 	ilink = ilink->next;
769421e10b8SBarry Smith       }
770421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
771421e10b8SBarry Smith     } else {
772421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
773421e10b8SBarry Smith       while (ilink) {
774421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
775421e10b8SBarry Smith 	ilink = ilink->next;
776421e10b8SBarry Smith       }
777421e10b8SBarry Smith     }
778421e10b8SBarry Smith   } else {
779421e10b8SBarry Smith     if (!jac->w1) {
780421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
781421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
782421e10b8SBarry Smith     }
783421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
784421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
785421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
786421e10b8SBarry Smith       while (ilink->next) {
787421e10b8SBarry Smith         ilink = ilink->next;
7889989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
789421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
790421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
791421e10b8SBarry Smith       }
792421e10b8SBarry Smith       while (ilink->previous) {
793421e10b8SBarry Smith         ilink = ilink->previous;
7949989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
795421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
796421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
797421e10b8SBarry Smith       }
798421e10b8SBarry Smith     } else {
799421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
800421e10b8SBarry Smith 	ilink = ilink->next;
801421e10b8SBarry Smith       }
802421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
803421e10b8SBarry Smith       while (ilink->previous) {
804421e10b8SBarry Smith 	ilink = ilink->previous;
8059989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
806421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
807421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
808421e10b8SBarry Smith       }
809421e10b8SBarry Smith     }
810421e10b8SBarry Smith   }
811421e10b8SBarry Smith   CHKMEMQ;
812421e10b8SBarry Smith   PetscFunctionReturn(0);
813421e10b8SBarry Smith }
814421e10b8SBarry Smith 
8150971522cSBarry Smith #undef __FUNCT__
816574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
817574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8180971522cSBarry Smith {
8190971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8200971522cSBarry Smith   PetscErrorCode    ierr;
8215a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8220971522cSBarry Smith 
8230971522cSBarry Smith   PetscFunctionBegin;
8245a9f2f41SSatish Balay   while (ilink) {
825574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
826fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
827fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
828fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
829fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
830b5787286SJed Brown     ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
8315a9f2f41SSatish Balay     next = ilink->next;
8325a9f2f41SSatish Balay     ilink = next;
8330971522cSBarry Smith   }
83405b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
835574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
836574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
837574deadeSBarry Smith   } else if (jac->mat) {
838574deadeSBarry Smith     jac->mat = PETSC_NULL;
839574deadeSBarry Smith   }
84097bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
84168dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8426bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8436bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8446bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8456bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8466bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8476bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8486bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
84963ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
850574deadeSBarry Smith   PetscFunctionReturn(0);
851574deadeSBarry Smith }
852574deadeSBarry Smith 
853574deadeSBarry Smith #undef __FUNCT__
854574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
855574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
856574deadeSBarry Smith {
857574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
858574deadeSBarry Smith   PetscErrorCode    ierr;
859574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
860574deadeSBarry Smith 
861574deadeSBarry Smith   PetscFunctionBegin;
862574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
863574deadeSBarry Smith   while (ilink) {
8646bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
865574deadeSBarry Smith     next = ilink->next;
866574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
867574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
8685d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
869574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
870574deadeSBarry Smith     ilink = next;
871574deadeSBarry Smith   }
872574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
873c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8747b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
8757b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
8767b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
8777b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
8787b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
8790971522cSBarry Smith   PetscFunctionReturn(0);
8800971522cSBarry Smith }
8810971522cSBarry Smith 
8820971522cSBarry Smith #undef __FUNCT__
8830971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8840971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8850971522cSBarry Smith {
8861b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8876c924f48SJed Brown   PetscInt        bs;
888bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8899dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8903b224e63SBarry Smith   PCCompositeType ctype;
8911b9fc7fcSBarry Smith 
8920971522cSBarry Smith   PetscFunctionBegin;
8931b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
894acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
89551f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
89651f519a2SBarry Smith   if (flg) {
89751f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
89851f519a2SBarry Smith   }
899704ba839SBarry Smith 
900435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
901c0adfefeSBarry Smith   if (stokes) {
902c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
903c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
904c0adfefeSBarry Smith   }
905c0adfefeSBarry Smith 
9063b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
9073b224e63SBarry Smith   if (flg) {
9083b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
9093b224e63SBarry Smith   }
910d32f9abdSBarry Smith 
911c30613efSMatthew Knepley   /* Only setup fields once */
912c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
913d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
914d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
9156c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9166c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
917d32f9abdSBarry Smith   }
918c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
919c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
920084e4875SJed 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);
921c5d2311dSJed Brown   }
9221b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9230971522cSBarry Smith   PetscFunctionReturn(0);
9240971522cSBarry Smith }
9250971522cSBarry Smith 
9260971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9270971522cSBarry Smith 
9280971522cSBarry Smith EXTERN_C_BEGIN
9290971522cSBarry Smith #undef __FUNCT__
9300971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9315d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
9320971522cSBarry Smith {
93397bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9340971522cSBarry Smith   PetscErrorCode    ierr;
9355a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
93669a612a9SBarry Smith   char              prefix[128];
9375d4c12cdSJungho Lee   PetscInt          i;
9380971522cSBarry Smith 
9390971522cSBarry Smith   PetscFunctionBegin;
9406c924f48SJed Brown   if (jac->splitdefined) {
9416c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9426c924f48SJed Brown     PetscFunctionReturn(0);
9436c924f48SJed Brown   }
94451f519a2SBarry Smith   for (i=0; i<n; i++) {
945e32f2f54SBarry 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);
946e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
94751f519a2SBarry Smith   }
948704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
949a04f6461SBarry Smith   if (splitname) {
950db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
951a04f6461SBarry Smith   } else {
952a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
953a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
954a04f6461SBarry Smith   }
955704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9565a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9575d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
9585d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
9595a9f2f41SSatish Balay   ilink->nfields = n;
9605a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9617adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9621cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9635a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9649005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
96569a612a9SBarry Smith 
966a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9675a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9680971522cSBarry Smith 
9690971522cSBarry Smith   if (!next) {
9705a9f2f41SSatish Balay     jac->head       = ilink;
97151f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9720971522cSBarry Smith   } else {
9730971522cSBarry Smith     while (next->next) {
9740971522cSBarry Smith       next = next->next;
9750971522cSBarry Smith     }
9765a9f2f41SSatish Balay     next->next      = ilink;
97751f519a2SBarry Smith     ilink->previous = next;
9780971522cSBarry Smith   }
9790971522cSBarry Smith   jac->nsplits++;
9800971522cSBarry Smith   PetscFunctionReturn(0);
9810971522cSBarry Smith }
9820971522cSBarry Smith EXTERN_C_END
9830971522cSBarry Smith 
984e69d4d44SBarry Smith EXTERN_C_BEGIN
985e69d4d44SBarry Smith #undef __FUNCT__
986e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9877087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
988e69d4d44SBarry Smith {
989e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
990e69d4d44SBarry Smith   PetscErrorCode ierr;
991e69d4d44SBarry Smith 
992e69d4d44SBarry Smith   PetscFunctionBegin;
993519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
994e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
995e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
99613e0d083SBarry Smith   if (n) *n = jac->nsplits;
997e69d4d44SBarry Smith   PetscFunctionReturn(0);
998e69d4d44SBarry Smith }
999e69d4d44SBarry Smith EXTERN_C_END
10000971522cSBarry Smith 
10010971522cSBarry Smith EXTERN_C_BEGIN
10020971522cSBarry Smith #undef __FUNCT__
100369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
10047087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
10050971522cSBarry Smith {
10060971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10070971522cSBarry Smith   PetscErrorCode    ierr;
10080971522cSBarry Smith   PetscInt          cnt = 0;
10095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
10100971522cSBarry Smith 
10110971522cSBarry Smith   PetscFunctionBegin;
10125d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
10135a9f2f41SSatish Balay   while (ilink) {
10145a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
10155a9f2f41SSatish Balay     ilink = ilink->next;
10160971522cSBarry Smith   }
10175d480477SMatthew 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);
101813e0d083SBarry Smith   if (n) *n = jac->nsplits;
10190971522cSBarry Smith   PetscFunctionReturn(0);
10200971522cSBarry Smith }
10210971522cSBarry Smith EXTERN_C_END
10220971522cSBarry Smith 
1023704ba839SBarry Smith EXTERN_C_BEGIN
1024704ba839SBarry Smith #undef __FUNCT__
1025704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10267087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1027704ba839SBarry Smith {
1028704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1029704ba839SBarry Smith   PetscErrorCode    ierr;
1030704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1031704ba839SBarry Smith   char              prefix[128];
1032704ba839SBarry Smith 
1033704ba839SBarry Smith   PetscFunctionBegin;
10346c924f48SJed Brown   if (jac->splitdefined) {
10356c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10366c924f48SJed Brown     PetscFunctionReturn(0);
10376c924f48SJed Brown   }
103816913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1039a04f6461SBarry Smith   if (splitname) {
1040db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1041a04f6461SBarry Smith   } else {
1042b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1043b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1044a04f6461SBarry Smith   }
10451ab39975SBarry Smith   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1046b5787286SJed Brown   ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
1047b5787286SJed Brown   ilink->is      = is;
1048b5787286SJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1049b5787286SJed Brown   ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1050b5787286SJed Brown   ilink->is_col  = is;
1051704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1052704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10531cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1054704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10559005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1056704ba839SBarry Smith 
1057a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1058704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1059704ba839SBarry Smith 
1060704ba839SBarry Smith   if (!next) {
1061704ba839SBarry Smith     jac->head       = ilink;
1062704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1063704ba839SBarry Smith   } else {
1064704ba839SBarry Smith     while (next->next) {
1065704ba839SBarry Smith       next = next->next;
1066704ba839SBarry Smith     }
1067704ba839SBarry Smith     next->next      = ilink;
1068704ba839SBarry Smith     ilink->previous = next;
1069704ba839SBarry Smith   }
1070704ba839SBarry Smith   jac->nsplits++;
1071704ba839SBarry Smith 
1072704ba839SBarry Smith   PetscFunctionReturn(0);
1073704ba839SBarry Smith }
1074704ba839SBarry Smith EXTERN_C_END
1075704ba839SBarry Smith 
10760971522cSBarry Smith #undef __FUNCT__
10770971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10780971522cSBarry Smith /*@
10790971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10800971522cSBarry Smith 
1081ad4df100SBarry Smith     Logically Collective on PC
10820971522cSBarry Smith 
10830971522cSBarry Smith     Input Parameters:
10840971522cSBarry Smith +   pc  - the preconditioner context
1085a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10860971522cSBarry Smith .   n - the number of fields in this split
1087db4c96c1SJed Brown -   fields - the fields in this split
10880971522cSBarry Smith 
10890971522cSBarry Smith     Level: intermediate
10900971522cSBarry Smith 
1091d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1092d32f9abdSBarry Smith 
1093d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1094d32f9abdSBarry 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
1095d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1096d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1097d32f9abdSBarry Smith 
1098db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1099db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1100db4c96c1SJed Brown 
11015d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
11025d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
11035d4c12cdSJungho Lee      available when this routine is called.
11045d4c12cdSJungho Lee 
1105d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
11060971522cSBarry Smith 
11070971522cSBarry Smith @*/
11085d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
11090971522cSBarry Smith {
11104ac538c5SBarry Smith   PetscErrorCode ierr;
11110971522cSBarry Smith 
11120971522cSBarry Smith   PetscFunctionBegin;
11130700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1114db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1115db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1116db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
11175d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
11180971522cSBarry Smith   PetscFunctionReturn(0);
11190971522cSBarry Smith }
11200971522cSBarry Smith 
11210971522cSBarry Smith #undef __FUNCT__
1122704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1123704ba839SBarry Smith /*@
1124704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1125704ba839SBarry Smith 
1126ad4df100SBarry Smith     Logically Collective on PC
1127704ba839SBarry Smith 
1128704ba839SBarry Smith     Input Parameters:
1129704ba839SBarry Smith +   pc  - the preconditioner context
1130a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1131db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1132704ba839SBarry Smith 
1133d32f9abdSBarry Smith 
1134a6ffb8dbSJed Brown     Notes:
1135a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1136a6ffb8dbSJed Brown 
1137db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1138db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1139d32f9abdSBarry Smith 
1140704ba839SBarry Smith     Level: intermediate
1141704ba839SBarry Smith 
1142704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1143704ba839SBarry Smith 
1144704ba839SBarry Smith @*/
11457087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1146704ba839SBarry Smith {
11474ac538c5SBarry Smith   PetscErrorCode ierr;
1148704ba839SBarry Smith 
1149704ba839SBarry Smith   PetscFunctionBegin;
11500700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11517b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1152db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11534ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1154704ba839SBarry Smith   PetscFunctionReturn(0);
1155704ba839SBarry Smith }
1156704ba839SBarry Smith 
1157704ba839SBarry Smith #undef __FUNCT__
115857a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
115957a9adfeSMatthew G Knepley /*@
116057a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
116157a9adfeSMatthew G Knepley 
116257a9adfeSMatthew G Knepley     Logically Collective on PC
116357a9adfeSMatthew G Knepley 
116457a9adfeSMatthew G Knepley     Input Parameters:
116557a9adfeSMatthew G Knepley +   pc  - the preconditioner context
116657a9adfeSMatthew G Knepley -   splitname - name of this split
116757a9adfeSMatthew G Knepley 
116857a9adfeSMatthew G Knepley     Output Parameter:
116957a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
117057a9adfeSMatthew G Knepley 
117157a9adfeSMatthew G Knepley     Level: intermediate
117257a9adfeSMatthew G Knepley 
117357a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
117457a9adfeSMatthew G Knepley 
117557a9adfeSMatthew G Knepley @*/
117657a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
117757a9adfeSMatthew G Knepley {
117857a9adfeSMatthew G Knepley   PetscErrorCode ierr;
117957a9adfeSMatthew G Knepley 
118057a9adfeSMatthew G Knepley   PetscFunctionBegin;
118157a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
118257a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
118357a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
118457a9adfeSMatthew G Knepley   {
118557a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
118657a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
118757a9adfeSMatthew G Knepley     PetscBool         found;
118857a9adfeSMatthew G Knepley 
118957a9adfeSMatthew G Knepley     *is = PETSC_NULL;
119057a9adfeSMatthew G Knepley     while(ilink) {
119157a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
119257a9adfeSMatthew G Knepley       if (found) {
119357a9adfeSMatthew G Knepley         *is = ilink->is;
119457a9adfeSMatthew G Knepley         break;
119557a9adfeSMatthew G Knepley       }
119657a9adfeSMatthew G Knepley       ilink = ilink->next;
119757a9adfeSMatthew G Knepley     }
119857a9adfeSMatthew G Knepley   }
119957a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
120057a9adfeSMatthew G Knepley }
120157a9adfeSMatthew G Knepley 
120257a9adfeSMatthew G Knepley #undef __FUNCT__
120351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
120451f519a2SBarry Smith /*@
120551f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
120651f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
120751f519a2SBarry Smith 
1208ad4df100SBarry Smith     Logically Collective on PC
120951f519a2SBarry Smith 
121051f519a2SBarry Smith     Input Parameters:
121151f519a2SBarry Smith +   pc  - the preconditioner context
121251f519a2SBarry Smith -   bs - the block size
121351f519a2SBarry Smith 
121451f519a2SBarry Smith     Level: intermediate
121551f519a2SBarry Smith 
121651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
121751f519a2SBarry Smith 
121851f519a2SBarry Smith @*/
12197087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
122051f519a2SBarry Smith {
12214ac538c5SBarry Smith   PetscErrorCode ierr;
122251f519a2SBarry Smith 
122351f519a2SBarry Smith   PetscFunctionBegin;
12240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1225c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
12264ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
122751f519a2SBarry Smith   PetscFunctionReturn(0);
122851f519a2SBarry Smith }
122951f519a2SBarry Smith 
123051f519a2SBarry Smith #undef __FUNCT__
123169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12320971522cSBarry Smith /*@C
123369a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12340971522cSBarry Smith 
123569a612a9SBarry Smith    Collective on KSP
12360971522cSBarry Smith 
12370971522cSBarry Smith    Input Parameter:
12380971522cSBarry Smith .  pc - the preconditioner context
12390971522cSBarry Smith 
12400971522cSBarry Smith    Output Parameters:
124113e0d083SBarry Smith +  n - the number of splits
124269a612a9SBarry Smith -  pc - the array of KSP contexts
12430971522cSBarry Smith 
12440971522cSBarry Smith    Note:
1245d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1246d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12470971522cSBarry Smith 
124869a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12490971522cSBarry Smith 
12500971522cSBarry Smith    Level: advanced
12510971522cSBarry Smith 
12520971522cSBarry Smith .seealso: PCFIELDSPLIT
12530971522cSBarry Smith @*/
12547087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12550971522cSBarry Smith {
12564ac538c5SBarry Smith   PetscErrorCode ierr;
12570971522cSBarry Smith 
12580971522cSBarry Smith   PetscFunctionBegin;
12590700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
126013e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12614ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12620971522cSBarry Smith   PetscFunctionReturn(0);
12630971522cSBarry Smith }
12640971522cSBarry Smith 
1265e69d4d44SBarry Smith #undef __FUNCT__
1266e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1267e69d4d44SBarry Smith /*@
1268e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1269a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1270e69d4d44SBarry Smith 
1271e69d4d44SBarry Smith     Collective on PC
1272e69d4d44SBarry Smith 
1273e69d4d44SBarry Smith     Input Parameters:
1274e69d4d44SBarry Smith +   pc  - the preconditioner context
1275fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1276084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1277084e4875SJed Brown 
1278e69d4d44SBarry Smith     Options Database:
1279084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1280e69d4d44SBarry Smith 
1281fd1303e9SJungho Lee     Notes:
1282fd1303e9SJungho Lee $    If ptype is
1283fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1284fd1303e9SJungho Lee $             to this function).
1285fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1286fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1287fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1288fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1289fd1303e9SJungho Lee $             preconditioner
1290fd1303e9SJungho Lee 
1291fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1292fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1293fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1294fd1303e9SJungho Lee 
1295fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1296fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1297fd1303e9SJungho Lee 
1298e69d4d44SBarry Smith     Level: intermediate
1299e69d4d44SBarry Smith 
1300fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1301e69d4d44SBarry Smith 
1302e69d4d44SBarry Smith @*/
13037087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1304e69d4d44SBarry Smith {
13054ac538c5SBarry Smith   PetscErrorCode ierr;
1306e69d4d44SBarry Smith 
1307e69d4d44SBarry Smith   PetscFunctionBegin;
13080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13094ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1310e69d4d44SBarry Smith   PetscFunctionReturn(0);
1311e69d4d44SBarry Smith }
1312e69d4d44SBarry Smith 
1313e69d4d44SBarry Smith EXTERN_C_BEGIN
1314e69d4d44SBarry Smith #undef __FUNCT__
1315e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
13167087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1317e69d4d44SBarry Smith {
1318e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1319084e4875SJed Brown   PetscErrorCode  ierr;
1320e69d4d44SBarry Smith 
1321e69d4d44SBarry Smith   PetscFunctionBegin;
1322084e4875SJed Brown   jac->schurpre = ptype;
1323084e4875SJed Brown   if (pre) {
13246bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1325084e4875SJed Brown     jac->schur_user = pre;
1326084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1327084e4875SJed Brown   }
1328e69d4d44SBarry Smith   PetscFunctionReturn(0);
1329e69d4d44SBarry Smith }
1330e69d4d44SBarry Smith EXTERN_C_END
1331e69d4d44SBarry Smith 
133230ad9308SMatthew Knepley #undef __FUNCT__
133330ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
133430ad9308SMatthew Knepley /*@C
13358c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
133630ad9308SMatthew Knepley 
133730ad9308SMatthew Knepley    Collective on KSP
133830ad9308SMatthew Knepley 
133930ad9308SMatthew Knepley    Input Parameter:
134030ad9308SMatthew Knepley .  pc - the preconditioner context
134130ad9308SMatthew Knepley 
134230ad9308SMatthew Knepley    Output Parameters:
1343a04f6461SBarry Smith +  A00 - the (0,0) block
1344a04f6461SBarry Smith .  A01 - the (0,1) block
1345a04f6461SBarry Smith .  A10 - the (1,0) block
1346a04f6461SBarry Smith -  A11 - the (1,1) block
134730ad9308SMatthew Knepley 
134830ad9308SMatthew Knepley    Level: advanced
134930ad9308SMatthew Knepley 
135030ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
135130ad9308SMatthew Knepley @*/
1352a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
135330ad9308SMatthew Knepley {
135430ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
135530ad9308SMatthew Knepley 
135630ad9308SMatthew Knepley   PetscFunctionBegin;
13570700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1358c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1359a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1360a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1361a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1362a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
136330ad9308SMatthew Knepley   PetscFunctionReturn(0);
136430ad9308SMatthew Knepley }
136530ad9308SMatthew Knepley 
136679416396SBarry Smith EXTERN_C_BEGIN
136779416396SBarry Smith #undef __FUNCT__
136879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
13697087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
137079416396SBarry Smith {
137179416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1372e69d4d44SBarry Smith   PetscErrorCode ierr;
137379416396SBarry Smith 
137479416396SBarry Smith   PetscFunctionBegin;
137579416396SBarry Smith   jac->type = type;
13763b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
13773b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
13783b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1379e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1380e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1381e69d4d44SBarry Smith 
13823b224e63SBarry Smith   } else {
13833b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
13843b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1385e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13869e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13873b224e63SBarry Smith   }
138879416396SBarry Smith   PetscFunctionReturn(0);
138979416396SBarry Smith }
139079416396SBarry Smith EXTERN_C_END
139179416396SBarry Smith 
139251f519a2SBarry Smith EXTERN_C_BEGIN
139351f519a2SBarry Smith #undef __FUNCT__
139451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13957087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
139651f519a2SBarry Smith {
139751f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
139851f519a2SBarry Smith 
139951f519a2SBarry Smith   PetscFunctionBegin;
140065e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
140165e19b50SBarry 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);
140251f519a2SBarry Smith   jac->bs = bs;
140351f519a2SBarry Smith   PetscFunctionReturn(0);
140451f519a2SBarry Smith }
140551f519a2SBarry Smith EXTERN_C_END
140651f519a2SBarry Smith 
140779416396SBarry Smith #undef __FUNCT__
140879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1409bc08b0f1SBarry Smith /*@
141079416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
141179416396SBarry Smith 
141279416396SBarry Smith    Collective on PC
141379416396SBarry Smith 
141479416396SBarry Smith    Input Parameter:
141579416396SBarry Smith .  pc - the preconditioner context
141681540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
141779416396SBarry Smith 
141879416396SBarry Smith    Options Database Key:
1419a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
142079416396SBarry Smith 
142179416396SBarry Smith    Level: Developer
142279416396SBarry Smith 
142379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
142479416396SBarry Smith 
142579416396SBarry Smith .seealso: PCCompositeSetType()
142679416396SBarry Smith 
142779416396SBarry Smith @*/
14287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
142979416396SBarry Smith {
14304ac538c5SBarry Smith   PetscErrorCode ierr;
143179416396SBarry Smith 
143279416396SBarry Smith   PetscFunctionBegin;
14330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14344ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
143579416396SBarry Smith   PetscFunctionReturn(0);
143679416396SBarry Smith }
143779416396SBarry Smith 
14380971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
14390971522cSBarry Smith /*MC
1440a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1441a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
14420971522cSBarry Smith 
1443edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1444edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
14450971522cSBarry Smith 
1446a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
144769a612a9SBarry Smith          and set the options directly on the resulting KSP object
14480971522cSBarry Smith 
14490971522cSBarry Smith    Level: intermediate
14500971522cSBarry Smith 
145179416396SBarry Smith    Options Database Keys:
145281540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
145381540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
145481540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
145581540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
145681540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1457e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1458435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
145979416396SBarry Smith 
14605d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
14615d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
14625d4c12cdSJungho Lee 
1463c8a0d604SMatthew G Knepley    Notes:
1464c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1465d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1466d32f9abdSBarry Smith 
1467d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1468d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1469d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1470d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1471d32f9abdSBarry Smith 
1472c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1473c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1474c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1475c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1476c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1477a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1478c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1479c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1480c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1481a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1482c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1483c8a0d604SMatthew G Knepley      factorization type is set using -pc_fieldsplit_schur_factorization_type <diag, lower, upper, full>. The full is shown above,
1484*5668aaf4SBarry Smith      diag gives
1485c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1486c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
1487*5668aaf4SBarry Smith      note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of
1488c8a0d604SMatthew G Knepley $              (  A00   0 )
1489c8a0d604SMatthew G Knepley $              (  A10   S )
1490c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1491c8a0d604SMatthew G Knepley $              ( A00 A01 )
1492c8a0d604SMatthew G Knepley $              (  0   S  )
1493c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1494e69d4d44SBarry Smith 
1495edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1496edf189efSBarry Smith      is used automatically for a second block.
1497edf189efSBarry Smith 
1498ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1499ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1500ff218e97SBarry Smith 
1501ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1502ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1503ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
15040716a85fSBarry Smith 
1505a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
15060971522cSBarry Smith 
15077e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1508e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
15090971522cSBarry Smith M*/
15100971522cSBarry Smith 
15110971522cSBarry Smith EXTERN_C_BEGIN
15120971522cSBarry Smith #undef __FUNCT__
15130971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
15147087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
15150971522cSBarry Smith {
15160971522cSBarry Smith   PetscErrorCode ierr;
15170971522cSBarry Smith   PC_FieldSplit  *jac;
15180971522cSBarry Smith 
15190971522cSBarry Smith   PetscFunctionBegin;
152038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
15210971522cSBarry Smith   jac->bs        = -1;
15220971522cSBarry Smith   jac->nsplits   = 0;
15233e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1524e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1525c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
152651f519a2SBarry Smith 
15270971522cSBarry Smith   pc->data     = (void*)jac;
15280971522cSBarry Smith 
15290971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1530421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
15310971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1532574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
15330971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
15340971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
15350971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
15360971522cSBarry Smith   pc->ops->applyrichardson   = 0;
15370971522cSBarry Smith 
153869a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
153969a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
15400971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
15410971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1542704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1543704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
154479416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
154579416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
154651f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
154751f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
15480971522cSBarry Smith   PetscFunctionReturn(0);
15490971522cSBarry Smith }
15500971522cSBarry Smith EXTERN_C_END
15510971522cSBarry Smith 
15520971522cSBarry Smith 
1553a541d17aSBarry Smith 
1554