xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision bdddcaaa462b6a2d5fbabba66a3b530d990c70d6)
1dba47a55SKris Buschelman 
2b45d2f2cSJed Brown #include <petsc-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};
6c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
7c5d2311dSJed Brown 
80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
90971522cSBarry Smith struct _PC_FieldSplitLink {
1069a612a9SBarry Smith   KSP               ksp;
110971522cSBarry Smith   Vec               x,y;
12db4c96c1SJed Brown   char              *splitname;
130971522cSBarry Smith   PetscInt          nfields;
145d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
151b9fc7fcSBarry Smith   VecScatter        sctx;
165d4c12cdSJungho Lee   IS                is,is_col;
1751f519a2SBarry Smith   PC_FieldSplitLink next,previous;
180971522cSBarry Smith };
190971522cSBarry Smith 
200971522cSBarry Smith typedef struct {
2168dd23aaSBarry Smith   PCCompositeType                    type;
22ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
23ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
24ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
2530ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
2630ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
2779416396SBarry Smith   Vec                                *x,*y,w1,w2;
28519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
29519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3030ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
31ace3abfcSBarry Smith   PetscBool                          issetup;
3230ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3330ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
3430ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
35a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
36084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
37084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
38c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
3930ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4097bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4163ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
42c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
430971522cSBarry Smith } PC_FieldSplit;
440971522cSBarry Smith 
4516913363SBarry Smith /*
4616913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
4716913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
4816913363SBarry Smith    PC you could change this.
4916913363SBarry Smith */
50084e4875SJed Brown 
51e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
52084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
53084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
54084e4875SJed Brown {
55084e4875SJed Brown   switch (jac->schurpre) {
56084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
57084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
58084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
59084e4875SJed Brown     default:
60084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
61084e4875SJed Brown   }
62084e4875SJed Brown }
63084e4875SJed Brown 
64084e4875SJed Brown 
650971522cSBarry Smith #undef __FUNCT__
660971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
670971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
680971522cSBarry Smith {
690971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
700971522cSBarry Smith   PetscErrorCode    ierr;
71ace3abfcSBarry Smith   PetscBool         iascii;
720971522cSBarry Smith   PetscInt          i,j;
735a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
740971522cSBarry Smith 
750971522cSBarry Smith   PetscFunctionBegin;
76251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
770971522cSBarry Smith   if (iascii) {
782c9966d7SBarry Smith     if (jac->bs > 0) {
7951f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
802c9966d7SBarry Smith     } else {
812c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
822c9966d7SBarry Smith     }
8369a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
840971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
850971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
861ab39975SBarry Smith       if (ilink->fields) {
870971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
8879416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
895a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9079416396SBarry Smith 	  if (j > 0) {
9179416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9279416396SBarry Smith 	  }
935a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
940971522cSBarry Smith 	}
950971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9679416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
971ab39975SBarry Smith       } else {
981ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
991ab39975SBarry Smith       }
1005a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1015a9f2f41SSatish Balay       ilink = ilink->next;
1020971522cSBarry Smith     }
1030971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1040971522cSBarry Smith   } else {
10565e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1060971522cSBarry Smith   }
1070971522cSBarry Smith   PetscFunctionReturn(0);
1080971522cSBarry Smith }
1090971522cSBarry Smith 
1100971522cSBarry Smith #undef __FUNCT__
1113b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1123b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1133b224e63SBarry Smith {
1143b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1153b224e63SBarry Smith   PetscErrorCode    ierr;
116ace3abfcSBarry Smith   PetscBool         iascii;
1173b224e63SBarry Smith   PetscInt          i,j;
1183b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1193b224e63SBarry Smith   KSP               ksp;
1203b224e63SBarry Smith 
1213b224e63SBarry Smith   PetscFunctionBegin;
122251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1233b224e63SBarry Smith   if (iascii) {
1242c9966d7SBarry Smith     if (jac->bs > 0) {
125c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1262c9966d7SBarry Smith     } else {
1272c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1282c9966d7SBarry Smith     }
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1303b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1313b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1323b224e63SBarry Smith       if (ilink->fields) {
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1343b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1353b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1363b224e63SBarry Smith 	  if (j > 0) {
1373b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1383b224e63SBarry Smith 	  }
1393b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1403b224e63SBarry Smith 	}
1413b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1423b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1433b224e63SBarry Smith       } else {
1443b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1453b224e63SBarry Smith       }
1463b224e63SBarry Smith       ilink = ilink->next;
1473b224e63SBarry Smith     }
148435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1493b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15012cae6f2SJed Brown     if (jac->schur) {
15112cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1523b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15312cae6f2SJed Brown     } else {
15412cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15512cae6f2SJed Brown     }
1563b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
157435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1583b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15912cae6f2SJed Brown     if (jac->kspschur) {
1603b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16112cae6f2SJed Brown     } else {
16212cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16312cae6f2SJed Brown     }
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1663b224e63SBarry Smith   } else {
16765e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1683b224e63SBarry Smith   }
1693b224e63SBarry Smith   PetscFunctionReturn(0);
1703b224e63SBarry Smith }
1713b224e63SBarry Smith 
1723b224e63SBarry Smith #undef __FUNCT__
1736c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1746c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1756c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1766c924f48SJed Brown {
1776c924f48SJed Brown   PetscErrorCode ierr;
1786c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1795d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
1805d4c12cdSJungho Lee   PetscBool      flg,flg_col;
1815d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
1826c924f48SJed Brown 
1836c924f48SJed Brown   PetscFunctionBegin;
1846c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1855d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
1866c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1876c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1886c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1895d4c12cdSJungho Lee     ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
1906c924f48SJed Brown     nfields = jac->bs;
19129499fbbSJungho Lee     nfields_col = jac->bs;
1926c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1935d4c12cdSJungho Lee     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
1946c924f48SJed Brown     if (!flg) break;
1955d4c12cdSJungho Lee     else if (flg && !flg_col) {
1965d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1975d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
1985d4c12cdSJungho Lee     }
1995d4c12cdSJungho Lee     else {
2005d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2015d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
2025d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
2035d4c12cdSJungho Lee     }
2046c924f48SJed Brown   }
2056c924f48SJed Brown   if (i > 0) {
2066c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
2076c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2086c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2096c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2106c924f48SJed Brown   }
2116c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2125d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2136c924f48SJed Brown   PetscFunctionReturn(0);
2146c924f48SJed Brown }
2156c924f48SJed Brown 
2166c924f48SJed Brown #undef __FUNCT__
21769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
21869a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2190971522cSBarry Smith {
2200971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2210971522cSBarry Smith   PetscErrorCode    ierr;
2225a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2236ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2246c924f48SJed Brown   PetscInt          i;
2250971522cSBarry Smith 
2260971522cSBarry Smith   PetscFunctionBegin;
227d32f9abdSBarry Smith   if (!ilink) {
228d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
229d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2307b62db95SJungho Lee       PetscInt     numFields, f;
2310784a22cSJed Brown       char         **fieldNames;
2327b62db95SJungho Lee       IS          *fields;
233e7c4fc90SDmitry Karpeev       DM          *dms;
234e7c4fc90SDmitry Karpeev       ierr = DMCreateDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);    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);
242e7c4fc90SDmitry Karpeev       if(dms) {
2438b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2447b62db95SJungho Lee         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2457b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
2467b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
247e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]); CHKERRQ(ierr);
2482fa5ba8aSJed Brown         }
2497b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
2508b8307b2SJed Brown       }
25166ffff09SJed Brown     } else {
252521d7252SBarry Smith       if (jac->bs <= 0) {
253704ba839SBarry Smith         if (pc->pmat) {
254521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
255704ba839SBarry Smith         } else {
256704ba839SBarry Smith           jac->bs = 1;
257704ba839SBarry Smith         }
258521d7252SBarry Smith       }
259d32f9abdSBarry Smith 
260acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2616ce1633cSBarry Smith       if (stokes) {
2626ce1633cSBarry Smith         IS       zerodiags,rest;
2636ce1633cSBarry Smith         PetscInt nmin,nmax;
2646ce1633cSBarry Smith 
2656ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2666ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2676ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2686ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2696ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
270fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
271fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2726ce1633cSBarry Smith       } else {
273d32f9abdSBarry Smith         if (!flg) {
274d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
275d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2766c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2776c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
278d32f9abdSBarry Smith         }
2796c924f48SJed Brown         if (flg || !jac->splitdefined) {
280d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
281db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2826c924f48SJed Brown             char splitname[8];
2836c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2845d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
28579416396SBarry Smith           }
2865d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
287521d7252SBarry Smith         }
28866ffff09SJed Brown       }
2896ce1633cSBarry Smith     }
290edf189efSBarry Smith   } else if (jac->nsplits == 1) {
291edf189efSBarry Smith     if (ilink->is) {
292edf189efSBarry Smith       IS       is2;
293edf189efSBarry Smith       PetscInt nmin,nmax;
294edf189efSBarry Smith 
295edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
296edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
297db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
298fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
2997b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
30063ec74ffSBarry Smith   } else if (jac->reset) {
30163ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
302d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
303d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
304d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
305d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
306d0af7cd3SBarry Smith       PetscBool dmcomposite;
307251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
308d0af7cd3SBarry Smith       if (dmcomposite) {
309d0af7cd3SBarry Smith         PetscInt nDM;
310d0af7cd3SBarry Smith         IS       *fields;
3117b62db95SJungho Lee         DM       *dms;
312d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
313d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
314d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3157b62db95SJungho Lee         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3167b62db95SJungho Lee         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
317d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3187b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3197b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
320d0af7cd3SBarry Smith           ilink->is = fields[i];
321d0af7cd3SBarry Smith           ilink     = ilink->next;
322edf189efSBarry Smith         }
323d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3247b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
325d0af7cd3SBarry Smith       }
326d0af7cd3SBarry Smith     } else {
327d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
328d0af7cd3SBarry Smith       if (stokes) {
329d0af7cd3SBarry Smith         IS       zerodiags,rest;
330d0af7cd3SBarry Smith         PetscInt nmin,nmax;
331d0af7cd3SBarry Smith 
332d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
333d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
334d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
33520b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
33620b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
337d0af7cd3SBarry Smith         ilink->is       = rest;
338d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
33920b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
340d0af7cd3SBarry Smith     }
341d0af7cd3SBarry Smith   }
342d0af7cd3SBarry Smith 
3437b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
34469a612a9SBarry Smith   PetscFunctionReturn(0);
34569a612a9SBarry Smith }
34669a612a9SBarry Smith 
34769a612a9SBarry Smith #undef __FUNCT__
34869a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
34969a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
35069a612a9SBarry Smith {
35169a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
35269a612a9SBarry Smith   PetscErrorCode    ierr;
3535a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
3542c9966d7SBarry Smith   PetscInt          i,nsplit;
35569a612a9SBarry Smith   MatStructure      flag = pc->flag;
3565f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
35769a612a9SBarry Smith 
35869a612a9SBarry Smith   PetscFunctionBegin;
35969a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
36097bbdb24SBarry Smith   nsplit = jac->nsplits;
3615a9f2f41SSatish Balay   ilink  = jac->head;
36297bbdb24SBarry Smith 
36397bbdb24SBarry Smith   /* get the matrices for each split */
364704ba839SBarry Smith   if (!jac->issetup) {
3651b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
36697bbdb24SBarry Smith 
367704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
368704ba839SBarry Smith 
3695d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
3702c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
3712c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
3722c9966d7SBarry Smith     }
37351f519a2SBarry Smith     bs     = jac->bs;
37497bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
3751b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3761b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3771b9fc7fcSBarry Smith       if (jac->defaultsplit) {
378704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
3795f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
380704ba839SBarry Smith       } else if (!ilink->is) {
381ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3825f4ab4e1SJungho Lee           PetscInt   *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
3835a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3845f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
3851b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3861b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3871b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
3885f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
38997bbdb24SBarry Smith             }
39097bbdb24SBarry Smith           }
391d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
3925f4ab4e1SJungho Lee           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
393ccb205f8SBarry Smith         } else {
394704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
3955f4ab4e1SJungho Lee           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
396ccb205f8SBarry Smith        }
3973e197d65SBarry Smith       }
398704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
3995f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
4005f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
401704ba839SBarry Smith       ilink = ilink->next;
4021b9fc7fcSBarry Smith     }
4031b9fc7fcSBarry Smith   }
4041b9fc7fcSBarry Smith 
405704ba839SBarry Smith   ilink  = jac->head;
40697bbdb24SBarry Smith   if (!jac->pmat) {
407*bdddcaaaSMatthew G Knepley     Vec xtmp;
408*bdddcaaaSMatthew G Knepley 
409*bdddcaaaSMatthew G Knepley     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
410cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
411*bdddcaaaSMatthew G Knepley     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
412cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
413*bdddcaaaSMatthew G Knepley       MatNullSpace sp;
414*bdddcaaaSMatthew G Knepley 
4155f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
416*bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
417*bdddcaaaSMatthew G Knepley       ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
418*bdddcaaaSMatthew G Knepley       ilink->x = jac->x[i]; ilink->y = jac->y[i];
419*bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
420*bdddcaaaSMatthew G Knepley       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
421*bdddcaaaSMatthew G Knepley       /* HACK: Check for the constant null space */
422*bdddcaaaSMatthew G Knepley       ierr = MatGetNullSpace(pc->pmat, &sp);CHKERRQ(ierr);
423*bdddcaaaSMatthew G Knepley       if (sp) {
424*bdddcaaaSMatthew G Knepley         MatNullSpace subsp;
425*bdddcaaaSMatthew G Knepley         Vec          ftmp, gtmp;
426*bdddcaaaSMatthew G Knepley         PetscScalar  norm;
427*bdddcaaaSMatthew G Knepley         PetscInt     N;
428*bdddcaaaSMatthew G Knepley 
429*bdddcaaaSMatthew G Knepley         ierr = MatGetVecs(pc->pmat,     &gtmp, PETSC_NULL);CHKERRQ(ierr);
430*bdddcaaaSMatthew G Knepley         ierr = MatGetVecs(jac->pmat[i], &ftmp, PETSC_NULL);CHKERRQ(ierr);
431*bdddcaaaSMatthew G Knepley         ierr = VecGetSize(ftmp, &N);CHKERRQ(ierr);
432*bdddcaaaSMatthew G Knepley         ierr = VecSet(ftmp, 1.0/N);CHKERRQ(ierr);
433*bdddcaaaSMatthew G Knepley         ierr = VecScatterBegin(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
434*bdddcaaaSMatthew G Knepley         ierr = VecScatterEnd(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
435*bdddcaaaSMatthew G Knepley         ierr = MatNullSpaceRemove(sp, gtmp, PETSC_NULL);CHKERRQ(ierr);
436*bdddcaaaSMatthew G Knepley         ierr = VecNorm(gtmp, NORM_2, &norm);CHKERRQ(ierr);
437*bdddcaaaSMatthew G Knepley         if (norm < 1.0e-10) {
438*bdddcaaaSMatthew G Knepley           ierr  = MatNullSpaceCreate(((PetscObject)pc)->comm, PETSC_TRUE, 0, PETSC_NULL, &subsp);CHKERRQ(ierr);
439*bdddcaaaSMatthew G Knepley           ierr  = MatSetNullSpace(jac->pmat[i], subsp);CHKERRQ(ierr);
440*bdddcaaaSMatthew G Knepley           ierr  = MatNullSpaceDestroy(&subsp);CHKERRQ(ierr);
441*bdddcaaaSMatthew G Knepley         }
442*bdddcaaaSMatthew G Knepley         ierr = VecDestroy(&ftmp);CHKERRQ(ierr);
443*bdddcaaaSMatthew G Knepley         ierr = VecDestroy(&gtmp);CHKERRQ(ierr);
444*bdddcaaaSMatthew G Knepley       }
445704ba839SBarry Smith       ilink = ilink->next;
446cf502942SBarry Smith     }
447*bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
44897bbdb24SBarry Smith   } else {
449cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4505f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
451704ba839SBarry Smith       ilink = ilink->next;
452cf502942SBarry Smith     }
45397bbdb24SBarry Smith   }
454519d70e2SJed Brown   if (jac->realdiagonal) {
455519d70e2SJed Brown     ilink = jac->head;
456519d70e2SJed Brown     if (!jac->mat) {
457519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
458519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4595f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
460519d70e2SJed Brown         ilink = ilink->next;
461519d70e2SJed Brown       }
462519d70e2SJed Brown     } else {
463519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4645f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
465519d70e2SJed Brown         ilink = ilink->next;
466519d70e2SJed Brown       }
467519d70e2SJed Brown     }
468519d70e2SJed Brown   } else {
469519d70e2SJed Brown     jac->mat = jac->pmat;
470519d70e2SJed Brown   }
47197bbdb24SBarry Smith 
4726c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
47368dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
47468dd23aaSBarry Smith     ilink  = jac->head;
47568dd23aaSBarry Smith     if (!jac->Afield) {
47668dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
47768dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4784aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
47968dd23aaSBarry Smith         ilink = ilink->next;
48068dd23aaSBarry Smith       }
48168dd23aaSBarry Smith     } else {
48268dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4834aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
48468dd23aaSBarry Smith         ilink = ilink->next;
48568dd23aaSBarry Smith       }
48668dd23aaSBarry Smith     }
48768dd23aaSBarry Smith   }
48868dd23aaSBarry Smith 
4893b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4903b224e63SBarry Smith     IS       ccis;
4914aa3045dSJed Brown     PetscInt rstart,rend;
492093c86ffSJed Brown     char     lscname[256];
493093c86ffSJed Brown     PetscObject LSC_L;
494e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
49568dd23aaSBarry Smith 
496e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
497e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
498e6cab6aaSJed Brown 
4993b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
5003b224e63SBarry Smith     if (jac->schur) {
5013b224e63SBarry Smith       ilink = jac->head;
50249bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5034aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
504fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5053b224e63SBarry Smith       ilink = ilink->next;
50649bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5074aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
508fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
509519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
510084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
5113b224e63SBarry Smith 
5123b224e63SBarry Smith      } else {
5131cee3971SBarry Smith       KSP ksp;
5146c924f48SJed Brown       char schurprefix[256];
515*bdddcaaaSMatthew G Knepley       MatNullSpace sp;
5163b224e63SBarry Smith 
517a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
5183b224e63SBarry Smith       ilink = jac->head;
51949bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5204aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
521fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5223b224e63SBarry Smith       ilink = ilink->next;
52349bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5244aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
525fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5267d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
527519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
528*bdddcaaaSMatthew G Knepley       ierr  = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr);
529*bdddcaaaSMatthew G Knepley       if (sp) {ierr  = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);}
530a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
5311cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
532e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
533a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
534a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
53520b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
53620b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
5377b62db95SJungho Lee       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
5383b224e63SBarry Smith 
5393b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
5409005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
5411cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
542084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
543084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
544084e4875SJed Brown         PC pc;
545cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
546084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
547084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
548e69d4d44SBarry Smith       }
5496c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5506c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5513b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
55220b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
55320b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5543b224e63SBarry Smith     }
555093c86ffSJed Brown 
556093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
557093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
558093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
559093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
560093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
561093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
562093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
563093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
564093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
5653b224e63SBarry Smith   } else {
56697bbdb24SBarry Smith     /* set up the individual PCs */
56797bbdb24SBarry Smith     i    = 0;
5685a9f2f41SSatish Balay     ilink = jac->head;
5695a9f2f41SSatish Balay     while (ilink) {
570519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5713b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
572c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5735a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
57497bbdb24SBarry Smith       i++;
5755a9f2f41SSatish Balay       ilink = ilink->next;
5760971522cSBarry Smith     }
5773b224e63SBarry Smith   }
5783b224e63SBarry Smith 
579c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5800971522cSBarry Smith   PetscFunctionReturn(0);
5810971522cSBarry Smith }
5820971522cSBarry Smith 
5835a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
584ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
585ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5865a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
587ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
588ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
58979416396SBarry Smith 
5900971522cSBarry Smith #undef __FUNCT__
5913b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5923b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5933b224e63SBarry Smith {
5943b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5953b224e63SBarry Smith   PetscErrorCode    ierr;
5963b224e63SBarry Smith   KSP               ksp;
5973b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5983b224e63SBarry Smith 
5993b224e63SBarry Smith   PetscFunctionBegin;
6003b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
6013b224e63SBarry Smith 
602c5d2311dSJed Brown   switch (jac->schurfactorization) {
603c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
604a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
605c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
606c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
607c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
608c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
609c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
610c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
611c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
615c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
616c5d2311dSJed Brown       break;
617c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
618a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
619c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
620c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
621c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
622c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
623c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
624c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
625c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
626c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
628c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
629c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
630c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
631c5d2311dSJed Brown       break;
632c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
633a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
634c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
635c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
636c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
637c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
638c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
639c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
640c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
641c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
642c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
643c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
644c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
645c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
646c5d2311dSJed Brown       break;
647c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_FULL:
648a04f6461SBarry 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 */
6493b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6503b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6513b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6523b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6533b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6543b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6553b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6563b224e63SBarry Smith 
6573b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6583b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6593b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6603b224e63SBarry Smith 
6613b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6623b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6633b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6643b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6653b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
666c5d2311dSJed Brown   }
6673b224e63SBarry Smith   PetscFunctionReturn(0);
6683b224e63SBarry Smith }
6693b224e63SBarry Smith 
6703b224e63SBarry Smith #undef __FUNCT__
6710971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6720971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6730971522cSBarry Smith {
6740971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6750971522cSBarry Smith   PetscErrorCode    ierr;
6765a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
677939b8a20SBarry Smith   PetscInt          cnt,bs;
6780971522cSBarry Smith 
6790971522cSBarry Smith   PetscFunctionBegin;
68051f519a2SBarry Smith   CHKMEMQ;
68179416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6821b9fc7fcSBarry Smith     if (jac->defaultsplit) {
683939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
684939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
685939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
686939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
6870971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6885a9f2f41SSatish Balay       while (ilink) {
6895a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6905a9f2f41SSatish Balay         ilink = ilink->next;
6910971522cSBarry Smith       }
6920971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6931b9fc7fcSBarry Smith     } else {
694efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6955a9f2f41SSatish Balay       while (ilink) {
6965a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6975a9f2f41SSatish Balay         ilink = ilink->next;
6981b9fc7fcSBarry Smith       }
6991b9fc7fcSBarry Smith     }
70016913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
70179416396SBarry Smith     if (!jac->w1) {
70279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
70379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
70479416396SBarry Smith     }
705efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
7065a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7073e197d65SBarry Smith     cnt = 1;
7085a9f2f41SSatish Balay     while (ilink->next) {
7095a9f2f41SSatish Balay       ilink = ilink->next;
7103e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7113e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7123e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7133e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7143e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7153e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7163e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7173e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7183e197d65SBarry Smith     }
71951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
72011755939SBarry Smith       cnt -= 2;
72151f519a2SBarry Smith       while (ilink->previous) {
72251f519a2SBarry Smith         ilink = ilink->previous;
72311755939SBarry Smith         /* compute the residual only over the part of the vector needed */
72411755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
72511755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
72611755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72711755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72811755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
72911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
73011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
73151f519a2SBarry Smith       }
73211755939SBarry Smith     }
73365e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
73451f519a2SBarry Smith   CHKMEMQ;
7350971522cSBarry Smith   PetscFunctionReturn(0);
7360971522cSBarry Smith }
7370971522cSBarry Smith 
738421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
739ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
740ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
741421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
742ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
743ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
744421e10b8SBarry Smith 
745421e10b8SBarry Smith #undef __FUNCT__
7468c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
747421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
748421e10b8SBarry Smith {
749421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
750421e10b8SBarry Smith   PetscErrorCode    ierr;
751421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
752939b8a20SBarry Smith   PetscInt          bs;
753421e10b8SBarry Smith 
754421e10b8SBarry Smith   PetscFunctionBegin;
755421e10b8SBarry Smith   CHKMEMQ;
756421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
757421e10b8SBarry Smith     if (jac->defaultsplit) {
758939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
759939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
760939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
761939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
762421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
763421e10b8SBarry Smith       while (ilink) {
764421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
765421e10b8SBarry Smith 	ilink = ilink->next;
766421e10b8SBarry Smith       }
767421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
768421e10b8SBarry Smith     } else {
769421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
770421e10b8SBarry Smith       while (ilink) {
771421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
772421e10b8SBarry Smith 	ilink = ilink->next;
773421e10b8SBarry Smith       }
774421e10b8SBarry Smith     }
775421e10b8SBarry Smith   } else {
776421e10b8SBarry Smith     if (!jac->w1) {
777421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
778421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
779421e10b8SBarry Smith     }
780421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
781421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
782421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
783421e10b8SBarry Smith       while (ilink->next) {
784421e10b8SBarry Smith         ilink = ilink->next;
7859989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
786421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
787421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
788421e10b8SBarry Smith       }
789421e10b8SBarry Smith       while (ilink->previous) {
790421e10b8SBarry Smith         ilink = ilink->previous;
7919989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
792421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
793421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
794421e10b8SBarry Smith       }
795421e10b8SBarry Smith     } else {
796421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
797421e10b8SBarry Smith 	ilink = ilink->next;
798421e10b8SBarry Smith       }
799421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
800421e10b8SBarry Smith       while (ilink->previous) {
801421e10b8SBarry Smith 	ilink = ilink->previous;
8029989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
803421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
804421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
805421e10b8SBarry Smith       }
806421e10b8SBarry Smith     }
807421e10b8SBarry Smith   }
808421e10b8SBarry Smith   CHKMEMQ;
809421e10b8SBarry Smith   PetscFunctionReturn(0);
810421e10b8SBarry Smith }
811421e10b8SBarry Smith 
8120971522cSBarry Smith #undef __FUNCT__
813574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
814574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8150971522cSBarry Smith {
8160971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8170971522cSBarry Smith   PetscErrorCode    ierr;
8185a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8190971522cSBarry Smith 
8200971522cSBarry Smith   PetscFunctionBegin;
8215a9f2f41SSatish Balay   while (ilink) {
822574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
823fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
824fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
825fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
826fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
827b5787286SJed Brown     ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
8285a9f2f41SSatish Balay     next = ilink->next;
8295a9f2f41SSatish Balay     ilink = next;
8300971522cSBarry Smith   }
83105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
832574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
833574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
834574deadeSBarry Smith   } else if (jac->mat) {
835574deadeSBarry Smith     jac->mat = PETSC_NULL;
836574deadeSBarry Smith   }
83797bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
83868dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8396bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8406bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8416bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8426bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8436bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8446bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8456bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
84663ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
847574deadeSBarry Smith   PetscFunctionReturn(0);
848574deadeSBarry Smith }
849574deadeSBarry Smith 
850574deadeSBarry Smith #undef __FUNCT__
851574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
852574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
853574deadeSBarry Smith {
854574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
855574deadeSBarry Smith   PetscErrorCode    ierr;
856574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
857574deadeSBarry Smith 
858574deadeSBarry Smith   PetscFunctionBegin;
859574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
860574deadeSBarry Smith   while (ilink) {
8616bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
862574deadeSBarry Smith     next = ilink->next;
863574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
864574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
8655d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
866574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
867574deadeSBarry Smith     ilink = next;
868574deadeSBarry Smith   }
869574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
870c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8717b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
8727b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
8737b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
8747b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
8757b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
876ab1df9f5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr);
877c9c6ffaaSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr);
8780971522cSBarry Smith   PetscFunctionReturn(0);
8790971522cSBarry Smith }
8800971522cSBarry Smith 
8810971522cSBarry Smith #undef __FUNCT__
8820971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8830971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8840971522cSBarry Smith {
8851b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8866c924f48SJed Brown   PetscInt        bs;
887bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8889dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8893b224e63SBarry Smith   PCCompositeType ctype;
890e7c4fc90SDmitry Karpeev   DM              ddm;
891e7c4fc90SDmitry Karpeev   char            ddm_name[1025];
8921b9fc7fcSBarry Smith 
8930971522cSBarry Smith   PetscFunctionBegin;
8941b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
895e7c4fc90SDmitry Karpeev   if(pc->dm) {
896e7c4fc90SDmitry Karpeev     /* Allow the user to request a decomposition DM by name */
897e7c4fc90SDmitry Karpeev     ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
898e7c4fc90SDmitry Karpeev     ierr = PetscOptionsString("-pc_fieldsplit_decomposition_dm", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
899e7c4fc90SDmitry Karpeev     if(flg) {
900e7c4fc90SDmitry Karpeev       ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
901e7c4fc90SDmitry Karpeev       if(!ddm) {
902e7c4fc90SDmitry Karpeev         SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
903e7c4fc90SDmitry Karpeev       }
904e7c4fc90SDmitry Karpeev       ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
905e7c4fc90SDmitry Karpeev       ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
906e7c4fc90SDmitry Karpeev     }
907e7c4fc90SDmitry Karpeev   }
908acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
90951f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
91051f519a2SBarry Smith   if (flg) {
91151f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
91251f519a2SBarry Smith   }
913704ba839SBarry Smith 
914435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
915c0adfefeSBarry Smith   if (stokes) {
916c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
917c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
918c0adfefeSBarry Smith   }
919c0adfefeSBarry Smith 
9203b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
9213b224e63SBarry Smith   if (flg) {
9223b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
9233b224e63SBarry Smith   }
924c30613efSMatthew Knepley   /* Only setup fields once */
925c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
926d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
927d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
9286c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9296c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
930d32f9abdSBarry Smith   }
931c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
932c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
933c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
934c9c6ffaaSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_fact_type","Which off-diagonal parts of the block factorization to use","PCFieldSplitSetSchurFactType",PCFieldSplitSchurFactTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
935084e4875SJed 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);
936c5d2311dSJed Brown   }
9371b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9380971522cSBarry Smith   PetscFunctionReturn(0);
9390971522cSBarry Smith }
9400971522cSBarry Smith 
9410971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9420971522cSBarry Smith 
9430971522cSBarry Smith EXTERN_C_BEGIN
9440971522cSBarry Smith #undef __FUNCT__
9450971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9465d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
9470971522cSBarry Smith {
94897bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9490971522cSBarry Smith   PetscErrorCode    ierr;
9505a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
95169a612a9SBarry Smith   char              prefix[128];
9525d4c12cdSJungho Lee   PetscInt          i;
9530971522cSBarry Smith 
9540971522cSBarry Smith   PetscFunctionBegin;
9556c924f48SJed Brown   if (jac->splitdefined) {
9566c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9576c924f48SJed Brown     PetscFunctionReturn(0);
9586c924f48SJed Brown   }
95951f519a2SBarry Smith   for (i=0; i<n; i++) {
960e32f2f54SBarry 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);
961e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
96251f519a2SBarry Smith   }
963704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
964a04f6461SBarry Smith   if (splitname) {
965db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
966a04f6461SBarry Smith   } else {
967a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
968a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
969a04f6461SBarry Smith   }
970704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9715a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9725d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
9735d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
9745a9f2f41SSatish Balay   ilink->nfields = n;
9755a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9767adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9771cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9785a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9799005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
98069a612a9SBarry Smith 
981a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9825a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9830971522cSBarry Smith 
9840971522cSBarry Smith   if (!next) {
9855a9f2f41SSatish Balay     jac->head       = ilink;
98651f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9870971522cSBarry Smith   } else {
9880971522cSBarry Smith     while (next->next) {
9890971522cSBarry Smith       next = next->next;
9900971522cSBarry Smith     }
9915a9f2f41SSatish Balay     next->next      = ilink;
99251f519a2SBarry Smith     ilink->previous = next;
9930971522cSBarry Smith   }
9940971522cSBarry Smith   jac->nsplits++;
9950971522cSBarry Smith   PetscFunctionReturn(0);
9960971522cSBarry Smith }
9970971522cSBarry Smith EXTERN_C_END
9980971522cSBarry Smith 
999e69d4d44SBarry Smith EXTERN_C_BEGIN
1000e69d4d44SBarry Smith #undef __FUNCT__
1001e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
10027087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1003e69d4d44SBarry Smith {
1004e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1005e69d4d44SBarry Smith   PetscErrorCode ierr;
1006e69d4d44SBarry Smith 
1007e69d4d44SBarry Smith   PetscFunctionBegin;
1008519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1009e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1010e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
101113e0d083SBarry Smith   if (n) *n = jac->nsplits;
1012e69d4d44SBarry Smith   PetscFunctionReturn(0);
1013e69d4d44SBarry Smith }
1014e69d4d44SBarry Smith EXTERN_C_END
10150971522cSBarry Smith 
10160971522cSBarry Smith EXTERN_C_BEGIN
10170971522cSBarry Smith #undef __FUNCT__
101869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
10197087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
10200971522cSBarry Smith {
10210971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10220971522cSBarry Smith   PetscErrorCode    ierr;
10230971522cSBarry Smith   PetscInt          cnt = 0;
10245a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
10250971522cSBarry Smith 
10260971522cSBarry Smith   PetscFunctionBegin;
10275d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
10285a9f2f41SSatish Balay   while (ilink) {
10295a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
10305a9f2f41SSatish Balay     ilink = ilink->next;
10310971522cSBarry Smith   }
10325d480477SMatthew 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);
103313e0d083SBarry Smith   if (n) *n = jac->nsplits;
10340971522cSBarry Smith   PetscFunctionReturn(0);
10350971522cSBarry Smith }
10360971522cSBarry Smith EXTERN_C_END
10370971522cSBarry Smith 
1038704ba839SBarry Smith EXTERN_C_BEGIN
1039704ba839SBarry Smith #undef __FUNCT__
1040704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10417087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1042704ba839SBarry Smith {
1043704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1044704ba839SBarry Smith   PetscErrorCode    ierr;
1045704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1046704ba839SBarry Smith   char              prefix[128];
1047704ba839SBarry Smith 
1048704ba839SBarry Smith   PetscFunctionBegin;
10496c924f48SJed Brown   if (jac->splitdefined) {
10506c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10516c924f48SJed Brown     PetscFunctionReturn(0);
10526c924f48SJed Brown   }
105316913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1054a04f6461SBarry Smith   if (splitname) {
1055db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1056a04f6461SBarry Smith   } else {
1057b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1058b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1059a04f6461SBarry Smith   }
10601ab39975SBarry Smith   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1061b5787286SJed Brown   ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
1062b5787286SJed Brown   ilink->is      = is;
1063b5787286SJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1064b5787286SJed Brown   ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1065b5787286SJed Brown   ilink->is_col  = is;
1066704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1067704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10681cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1069704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10709005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1071704ba839SBarry Smith 
1072a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1073704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1074704ba839SBarry Smith 
1075704ba839SBarry Smith   if (!next) {
1076704ba839SBarry Smith     jac->head       = ilink;
1077704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1078704ba839SBarry Smith   } else {
1079704ba839SBarry Smith     while (next->next) {
1080704ba839SBarry Smith       next = next->next;
1081704ba839SBarry Smith     }
1082704ba839SBarry Smith     next->next      = ilink;
1083704ba839SBarry Smith     ilink->previous = next;
1084704ba839SBarry Smith   }
1085704ba839SBarry Smith   jac->nsplits++;
1086704ba839SBarry Smith 
1087704ba839SBarry Smith   PetscFunctionReturn(0);
1088704ba839SBarry Smith }
1089704ba839SBarry Smith EXTERN_C_END
1090704ba839SBarry Smith 
10910971522cSBarry Smith #undef __FUNCT__
10920971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10930971522cSBarry Smith /*@
10940971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10950971522cSBarry Smith 
1096ad4df100SBarry Smith     Logically Collective on PC
10970971522cSBarry Smith 
10980971522cSBarry Smith     Input Parameters:
10990971522cSBarry Smith +   pc  - the preconditioner context
1100a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
11010971522cSBarry Smith .   n - the number of fields in this split
1102db4c96c1SJed Brown -   fields - the fields in this split
11030971522cSBarry Smith 
11040971522cSBarry Smith     Level: intermediate
11050971522cSBarry Smith 
1106d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1107d32f9abdSBarry Smith 
1108d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1109d32f9abdSBarry 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
1110d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1111d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1112d32f9abdSBarry Smith 
1113db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1114db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1115db4c96c1SJed Brown 
11165d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
11175d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
11185d4c12cdSJungho Lee      available when this routine is called.
11195d4c12cdSJungho Lee 
1120d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
11210971522cSBarry Smith 
11220971522cSBarry Smith @*/
11235d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
11240971522cSBarry Smith {
11254ac538c5SBarry Smith   PetscErrorCode ierr;
11260971522cSBarry Smith 
11270971522cSBarry Smith   PetscFunctionBegin;
11280700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1129db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1130db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1131db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
11325d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
11330971522cSBarry Smith   PetscFunctionReturn(0);
11340971522cSBarry Smith }
11350971522cSBarry Smith 
11360971522cSBarry Smith #undef __FUNCT__
1137704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1138704ba839SBarry Smith /*@
1139704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1140704ba839SBarry Smith 
1141ad4df100SBarry Smith     Logically Collective on PC
1142704ba839SBarry Smith 
1143704ba839SBarry Smith     Input Parameters:
1144704ba839SBarry Smith +   pc  - the preconditioner context
1145a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1146db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1147704ba839SBarry Smith 
1148d32f9abdSBarry Smith 
1149a6ffb8dbSJed Brown     Notes:
1150a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1151a6ffb8dbSJed Brown 
1152db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1153db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1154d32f9abdSBarry Smith 
1155704ba839SBarry Smith     Level: intermediate
1156704ba839SBarry Smith 
1157704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1158704ba839SBarry Smith 
1159704ba839SBarry Smith @*/
11607087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1161704ba839SBarry Smith {
11624ac538c5SBarry Smith   PetscErrorCode ierr;
1163704ba839SBarry Smith 
1164704ba839SBarry Smith   PetscFunctionBegin;
11650700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11667b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1167db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11684ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1169704ba839SBarry Smith   PetscFunctionReturn(0);
1170704ba839SBarry Smith }
1171704ba839SBarry Smith 
1172704ba839SBarry Smith #undef __FUNCT__
117357a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
117457a9adfeSMatthew G Knepley /*@
117557a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
117657a9adfeSMatthew G Knepley 
117757a9adfeSMatthew G Knepley     Logically Collective on PC
117857a9adfeSMatthew G Knepley 
117957a9adfeSMatthew G Knepley     Input Parameters:
118057a9adfeSMatthew G Knepley +   pc  - the preconditioner context
118157a9adfeSMatthew G Knepley -   splitname - name of this split
118257a9adfeSMatthew G Knepley 
118357a9adfeSMatthew G Knepley     Output Parameter:
118457a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
118557a9adfeSMatthew G Knepley 
118657a9adfeSMatthew G Knepley     Level: intermediate
118757a9adfeSMatthew G Knepley 
118857a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
118957a9adfeSMatthew G Knepley 
119057a9adfeSMatthew G Knepley @*/
119157a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
119257a9adfeSMatthew G Knepley {
119357a9adfeSMatthew G Knepley   PetscErrorCode ierr;
119457a9adfeSMatthew G Knepley 
119557a9adfeSMatthew G Knepley   PetscFunctionBegin;
119657a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119757a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
119857a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
119957a9adfeSMatthew G Knepley   {
120057a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
120157a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
120257a9adfeSMatthew G Knepley     PetscBool         found;
120357a9adfeSMatthew G Knepley 
120457a9adfeSMatthew G Knepley     *is = PETSC_NULL;
120557a9adfeSMatthew G Knepley     while(ilink) {
120657a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
120757a9adfeSMatthew G Knepley       if (found) {
120857a9adfeSMatthew G Knepley         *is = ilink->is;
120957a9adfeSMatthew G Knepley         break;
121057a9adfeSMatthew G Knepley       }
121157a9adfeSMatthew G Knepley       ilink = ilink->next;
121257a9adfeSMatthew G Knepley     }
121357a9adfeSMatthew G Knepley   }
121457a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
121557a9adfeSMatthew G Knepley }
121657a9adfeSMatthew G Knepley 
121757a9adfeSMatthew G Knepley #undef __FUNCT__
121851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
121951f519a2SBarry Smith /*@
122051f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
122151f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
122251f519a2SBarry Smith 
1223ad4df100SBarry Smith     Logically Collective on PC
122451f519a2SBarry Smith 
122551f519a2SBarry Smith     Input Parameters:
122651f519a2SBarry Smith +   pc  - the preconditioner context
122751f519a2SBarry Smith -   bs - the block size
122851f519a2SBarry Smith 
122951f519a2SBarry Smith     Level: intermediate
123051f519a2SBarry Smith 
123151f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
123251f519a2SBarry Smith 
123351f519a2SBarry Smith @*/
12347087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
123551f519a2SBarry Smith {
12364ac538c5SBarry Smith   PetscErrorCode ierr;
123751f519a2SBarry Smith 
123851f519a2SBarry Smith   PetscFunctionBegin;
12390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1240c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
12414ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
124251f519a2SBarry Smith   PetscFunctionReturn(0);
124351f519a2SBarry Smith }
124451f519a2SBarry Smith 
124551f519a2SBarry Smith #undef __FUNCT__
124669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12470971522cSBarry Smith /*@C
124869a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12490971522cSBarry Smith 
125069a612a9SBarry Smith    Collective on KSP
12510971522cSBarry Smith 
12520971522cSBarry Smith    Input Parameter:
12530971522cSBarry Smith .  pc - the preconditioner context
12540971522cSBarry Smith 
12550971522cSBarry Smith    Output Parameters:
125613e0d083SBarry Smith +  n - the number of splits
125769a612a9SBarry Smith -  pc - the array of KSP contexts
12580971522cSBarry Smith 
12590971522cSBarry Smith    Note:
1260d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1261d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12620971522cSBarry Smith 
126369a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12640971522cSBarry Smith 
12650971522cSBarry Smith    Level: advanced
12660971522cSBarry Smith 
12670971522cSBarry Smith .seealso: PCFIELDSPLIT
12680971522cSBarry Smith @*/
12697087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12700971522cSBarry Smith {
12714ac538c5SBarry Smith   PetscErrorCode ierr;
12720971522cSBarry Smith 
12730971522cSBarry Smith   PetscFunctionBegin;
12740700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
127513e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12764ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12770971522cSBarry Smith   PetscFunctionReturn(0);
12780971522cSBarry Smith }
12790971522cSBarry Smith 
1280e69d4d44SBarry Smith #undef __FUNCT__
1281e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1282e69d4d44SBarry Smith /*@
1283e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1284a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1285e69d4d44SBarry Smith 
1286e69d4d44SBarry Smith     Collective on PC
1287e69d4d44SBarry Smith 
1288e69d4d44SBarry Smith     Input Parameters:
1289e69d4d44SBarry Smith +   pc  - the preconditioner context
1290fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1291084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1292084e4875SJed Brown 
1293e69d4d44SBarry Smith     Options Database:
1294084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1295e69d4d44SBarry Smith 
1296fd1303e9SJungho Lee     Notes:
1297fd1303e9SJungho Lee $    If ptype is
1298fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1299fd1303e9SJungho Lee $             to this function).
1300fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1301fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1302fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1303fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1304fd1303e9SJungho Lee $             preconditioner
1305fd1303e9SJungho Lee 
1306fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1307fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1308fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1309fd1303e9SJungho Lee 
1310fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1311fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1312fd1303e9SJungho Lee 
1313e69d4d44SBarry Smith     Level: intermediate
1314e69d4d44SBarry Smith 
1315fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1316e69d4d44SBarry Smith 
1317e69d4d44SBarry Smith @*/
13187087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1319e69d4d44SBarry Smith {
13204ac538c5SBarry Smith   PetscErrorCode ierr;
1321e69d4d44SBarry Smith 
1322e69d4d44SBarry Smith   PetscFunctionBegin;
13230700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13244ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1325e69d4d44SBarry Smith   PetscFunctionReturn(0);
1326e69d4d44SBarry Smith }
1327e69d4d44SBarry Smith 
1328e69d4d44SBarry Smith EXTERN_C_BEGIN
1329e69d4d44SBarry Smith #undef __FUNCT__
1330e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
13317087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1332e69d4d44SBarry Smith {
1333e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1334084e4875SJed Brown   PetscErrorCode  ierr;
1335e69d4d44SBarry Smith 
1336e69d4d44SBarry Smith   PetscFunctionBegin;
1337084e4875SJed Brown   jac->schurpre = ptype;
1338084e4875SJed Brown   if (pre) {
13396bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1340084e4875SJed Brown     jac->schur_user = pre;
1341084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1342084e4875SJed Brown   }
1343e69d4d44SBarry Smith   PetscFunctionReturn(0);
1344e69d4d44SBarry Smith }
1345e69d4d44SBarry Smith EXTERN_C_END
1346e69d4d44SBarry Smith 
134730ad9308SMatthew Knepley #undef __FUNCT__
1348c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1349ab1df9f5SJed Brown /*@
1350c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1351ab1df9f5SJed Brown 
1352ab1df9f5SJed Brown     Collective on PC
1353ab1df9f5SJed Brown 
1354ab1df9f5SJed Brown     Input Parameters:
1355ab1df9f5SJed Brown +   pc  - the preconditioner context
1356c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1357ab1df9f5SJed Brown 
1358ab1df9f5SJed Brown     Options Database:
1359c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1360ab1df9f5SJed Brown 
1361ab1df9f5SJed Brown 
1362ab1df9f5SJed Brown     Level: intermediate
1363ab1df9f5SJed Brown 
1364ab1df9f5SJed Brown     Notes:
1365ab1df9f5SJed Brown     The FULL factorization is
1366ab1df9f5SJed Brown 
1367ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1368ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1369ab1df9f5SJed Brown 
13706be4592eSBarry Smith     where S = D - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D,
13716be4592eSBarry Smith     and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES).
1372ab1df9f5SJed Brown 
13736be4592eSBarry Smith     If applied exactly, FULL factorization is a direct solver. The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial
13746be4592eSBarry Smith     of degree 2, so KSPGMRES converges in 2 iterations. If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
13756be4592eSBarry Smith     application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice. With DIAG,
13766be4592eSBarry Smith     the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.
1377ab1df9f5SJed Brown 
13786be4592eSBarry Smith     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES. Note that a flexible method like KSPFGMRES
13796be4592eSBarry Smith     or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used inside a split).
1380ab1df9f5SJed Brown 
1381ab1df9f5SJed Brown     References:
1382ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1383ab1df9f5SJed Brown 
1384ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1385ab1df9f5SJed Brown 
1386ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1387ab1df9f5SJed Brown @*/
1388c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1389ab1df9f5SJed Brown {
1390ab1df9f5SJed Brown   PetscErrorCode ierr;
1391ab1df9f5SJed Brown 
1392ab1df9f5SJed Brown   PetscFunctionBegin;
1393ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1394c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1395ab1df9f5SJed Brown   PetscFunctionReturn(0);
1396ab1df9f5SJed Brown }
1397ab1df9f5SJed Brown 
1398ab1df9f5SJed Brown #undef __FUNCT__
1399c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1400c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1401ab1df9f5SJed Brown {
1402ab1df9f5SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1403ab1df9f5SJed Brown 
1404ab1df9f5SJed Brown   PetscFunctionBegin;
1405ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1406ab1df9f5SJed Brown   PetscFunctionReturn(0);
1407ab1df9f5SJed Brown }
1408ab1df9f5SJed Brown 
1409ab1df9f5SJed Brown #undef __FUNCT__
141030ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
141130ad9308SMatthew Knepley /*@C
14128c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
141330ad9308SMatthew Knepley 
141430ad9308SMatthew Knepley    Collective on KSP
141530ad9308SMatthew Knepley 
141630ad9308SMatthew Knepley    Input Parameter:
141730ad9308SMatthew Knepley .  pc - the preconditioner context
141830ad9308SMatthew Knepley 
141930ad9308SMatthew Knepley    Output Parameters:
1420a04f6461SBarry Smith +  A00 - the (0,0) block
1421a04f6461SBarry Smith .  A01 - the (0,1) block
1422a04f6461SBarry Smith .  A10 - the (1,0) block
1423a04f6461SBarry Smith -  A11 - the (1,1) block
142430ad9308SMatthew Knepley 
142530ad9308SMatthew Knepley    Level: advanced
142630ad9308SMatthew Knepley 
142730ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
142830ad9308SMatthew Knepley @*/
1429a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
143030ad9308SMatthew Knepley {
143130ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
143230ad9308SMatthew Knepley 
143330ad9308SMatthew Knepley   PetscFunctionBegin;
14340700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1435c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1436a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1437a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1438a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1439a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
144030ad9308SMatthew Knepley   PetscFunctionReturn(0);
144130ad9308SMatthew Knepley }
144230ad9308SMatthew Knepley 
144379416396SBarry Smith EXTERN_C_BEGIN
144479416396SBarry Smith #undef __FUNCT__
144579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
14467087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
144779416396SBarry Smith {
144879416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1449e69d4d44SBarry Smith   PetscErrorCode ierr;
145079416396SBarry Smith 
145179416396SBarry Smith   PetscFunctionBegin;
145279416396SBarry Smith   jac->type = type;
14533b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
14543b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
14553b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1456e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1457e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1458c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1459e69d4d44SBarry Smith 
14603b224e63SBarry Smith   } else {
14613b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
14623b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1463e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14649e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1465c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
14663b224e63SBarry Smith   }
146779416396SBarry Smith   PetscFunctionReturn(0);
146879416396SBarry Smith }
146979416396SBarry Smith EXTERN_C_END
147079416396SBarry Smith 
147151f519a2SBarry Smith EXTERN_C_BEGIN
147251f519a2SBarry Smith #undef __FUNCT__
147351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
14747087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
147551f519a2SBarry Smith {
147651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
147751f519a2SBarry Smith 
147851f519a2SBarry Smith   PetscFunctionBegin;
147965e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
148065e19b50SBarry 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);
148151f519a2SBarry Smith   jac->bs = bs;
148251f519a2SBarry Smith   PetscFunctionReturn(0);
148351f519a2SBarry Smith }
148451f519a2SBarry Smith EXTERN_C_END
148551f519a2SBarry Smith 
148679416396SBarry Smith #undef __FUNCT__
148779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1488bc08b0f1SBarry Smith /*@
148979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
149079416396SBarry Smith 
149179416396SBarry Smith    Collective on PC
149279416396SBarry Smith 
149379416396SBarry Smith    Input Parameter:
149479416396SBarry Smith .  pc - the preconditioner context
149581540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
149679416396SBarry Smith 
149779416396SBarry Smith    Options Database Key:
1498a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
149979416396SBarry Smith 
150079416396SBarry Smith    Level: Developer
150179416396SBarry Smith 
150279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
150379416396SBarry Smith 
150479416396SBarry Smith .seealso: PCCompositeSetType()
150579416396SBarry Smith 
150679416396SBarry Smith @*/
15077087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
150879416396SBarry Smith {
15094ac538c5SBarry Smith   PetscErrorCode ierr;
151079416396SBarry Smith 
151179416396SBarry Smith   PetscFunctionBegin;
15120700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15134ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
151479416396SBarry Smith   PetscFunctionReturn(0);
151579416396SBarry Smith }
151679416396SBarry Smith 
15170971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
15180971522cSBarry Smith /*MC
1519a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1520a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
15210971522cSBarry Smith 
1522edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1523edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
15240971522cSBarry Smith 
1525a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
152669a612a9SBarry Smith          and set the options directly on the resulting KSP object
15270971522cSBarry Smith 
15280971522cSBarry Smith    Level: intermediate
15290971522cSBarry Smith 
153079416396SBarry Smith    Options Database Keys:
153181540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
153281540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
153381540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
153481540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
153581540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1536e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1537435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
153879416396SBarry Smith 
15395d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
15405d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
15415d4c12cdSJungho Lee 
1542c8a0d604SMatthew G Knepley    Notes:
1543c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1544d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1545d32f9abdSBarry Smith 
1546d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1547d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1548d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1549d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1550d32f9abdSBarry Smith 
1551c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1552c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1553c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1554c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1555c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1556a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1557c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1558c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1559c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1560a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1561c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1562c9c6ffaaSJed Brown      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
15635668aaf4SBarry Smith      diag gives
1564c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1565c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
15665668aaf4SBarry 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
1567c8a0d604SMatthew G Knepley $              (  A00   0 )
1568c8a0d604SMatthew G Knepley $              (  A10   S )
1569c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1570c8a0d604SMatthew G Knepley $              ( A00 A01 )
1571c8a0d604SMatthew G Knepley $              (  0   S  )
1572c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1573e69d4d44SBarry Smith 
1574edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1575edf189efSBarry Smith      is used automatically for a second block.
1576edf189efSBarry Smith 
1577ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1578ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1579ff218e97SBarry Smith 
1580ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1581ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1582ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
15830716a85fSBarry Smith 
1584a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
15850971522cSBarry Smith 
15867e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1587e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
15880971522cSBarry Smith M*/
15890971522cSBarry Smith 
15900971522cSBarry Smith EXTERN_C_BEGIN
15910971522cSBarry Smith #undef __FUNCT__
15920971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
15937087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
15940971522cSBarry Smith {
15950971522cSBarry Smith   PetscErrorCode ierr;
15960971522cSBarry Smith   PC_FieldSplit  *jac;
15970971522cSBarry Smith 
15980971522cSBarry Smith   PetscFunctionBegin;
159938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
16000971522cSBarry Smith   jac->bs        = -1;
16010971522cSBarry Smith   jac->nsplits   = 0;
16023e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1603e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1604c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
160551f519a2SBarry Smith 
16060971522cSBarry Smith   pc->data     = (void*)jac;
16070971522cSBarry Smith 
16080971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1609421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
16100971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1611574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
16120971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
16130971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
16140971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
16150971522cSBarry Smith   pc->ops->applyrichardson   = 0;
16160971522cSBarry Smith 
161769a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
161869a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
16190971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
16200971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1621704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1622704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
162379416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
162479416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
162551f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
162651f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
16270971522cSBarry Smith   PetscFunctionReturn(0);
16280971522cSBarry Smith }
16290971522cSBarry Smith EXTERN_C_END
16300971522cSBarry Smith 
16310971522cSBarry Smith 
1632a541d17aSBarry Smith 
1633