xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 8caf3d721b8a28a7e6d895d5ad47be9517ebf833)
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     }
83a3df900dSMatthew G Knepley     if (jac->realdiagonal) {
84a3df900dSMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr);
85a3df900dSMatthew G Knepley     }
8669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
870971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
880971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
891ab39975SBarry Smith       if (ilink->fields) {
900971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9179416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
925a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9379416396SBarry Smith 	  if (j > 0) {
9479416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9579416396SBarry Smith 	  }
965a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
970971522cSBarry Smith 	}
980971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1001ab39975SBarry Smith       } else {
1011ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1021ab39975SBarry Smith       }
1035a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1045a9f2f41SSatish Balay       ilink = ilink->next;
1050971522cSBarry Smith     }
1060971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1070971522cSBarry Smith   } else {
10865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1090971522cSBarry Smith   }
1100971522cSBarry Smith   PetscFunctionReturn(0);
1110971522cSBarry Smith }
1120971522cSBarry Smith 
1130971522cSBarry Smith #undef __FUNCT__
1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1163b224e63SBarry Smith {
1173b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1183b224e63SBarry Smith   PetscErrorCode    ierr;
119ace3abfcSBarry Smith   PetscBool         iascii;
1203b224e63SBarry Smith   PetscInt          i,j;
1213b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1223b224e63SBarry Smith   KSP               ksp;
1233b224e63SBarry Smith 
1243b224e63SBarry Smith   PetscFunctionBegin;
125251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1263b224e63SBarry Smith   if (iascii) {
1272c9966d7SBarry Smith     if (jac->bs > 0) {
128c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1292c9966d7SBarry Smith     } else {
1302c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1312c9966d7SBarry Smith     }
132a3df900dSMatthew G Knepley     if (jac->realdiagonal) {
133a3df900dSMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr);
134a3df900dSMatthew G Knepley     }
1353e8b8b31SMatthew G Knepley     switch(jac->schurpre) {
1363e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
1373e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break;
138b02e2d75SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_DIAG:
1393e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break;
1403e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1413e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1423e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1433e8b8b31SMatthew G Knepley       } else {
1443e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);
1453e8b8b31SMatthew G Knepley       }
1463e8b8b31SMatthew G Knepley       break;
1473e8b8b31SMatthew G Knepley     default:
1483e8b8b31SMatthew G Knepley       SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1493e8b8b31SMatthew G Knepley     }
1503b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1513b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1523b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1533b224e63SBarry Smith       if (ilink->fields) {
1543b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1553b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1563b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1573b224e63SBarry Smith 	  if (j > 0) {
1583b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1593b224e63SBarry Smith 	  }
1603b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1613b224e63SBarry Smith 	}
1623b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1633b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1643b224e63SBarry Smith       } else {
1653b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1663b224e63SBarry Smith       }
1673b224e63SBarry Smith       ilink = ilink->next;
1683b224e63SBarry Smith     }
169435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1703b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17112cae6f2SJed Brown     if (jac->schur) {
17212cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1733b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
17412cae6f2SJed Brown     } else {
17512cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
17612cae6f2SJed Brown     }
1773b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
178435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1793b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18012cae6f2SJed Brown     if (jac->kspschur) {
1813b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
18212cae6f2SJed Brown     } else {
18312cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
18412cae6f2SJed Brown     }
1853b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1863b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1873b224e63SBarry Smith   } else {
18865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1893b224e63SBarry Smith   }
1903b224e63SBarry Smith   PetscFunctionReturn(0);
1913b224e63SBarry Smith }
1923b224e63SBarry Smith 
1933b224e63SBarry Smith #undef __FUNCT__
1946c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1956c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1966c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1976c924f48SJed Brown {
1986c924f48SJed Brown   PetscErrorCode ierr;
1996c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2005d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
2015d4c12cdSJungho Lee   PetscBool      flg,flg_col;
2025d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
2036c924f48SJed Brown 
2046c924f48SJed Brown   PetscFunctionBegin;
2056c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
2065d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
2076c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
208*8caf3d72SBarry Smith     ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
209*8caf3d72SBarry Smith     ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
210*8caf3d72SBarry Smith     ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
2116c924f48SJed Brown     nfields = jac->bs;
21229499fbbSJungho Lee     nfields_col = jac->bs;
2136c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
2145d4c12cdSJungho Lee     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
2156c924f48SJed Brown     if (!flg) break;
2165d4c12cdSJungho Lee     else if (flg && !flg_col) {
2175d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2185d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
2195d4c12cdSJungho Lee     }
2205d4c12cdSJungho Lee     else {
2215d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2225d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
2235d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
2245d4c12cdSJungho Lee     }
2256c924f48SJed Brown   }
2266c924f48SJed Brown   if (i > 0) {
2276c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
2286c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2296c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2306c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2316c924f48SJed Brown   }
2326c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2335d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2346c924f48SJed Brown   PetscFunctionReturn(0);
2356c924f48SJed Brown }
2366c924f48SJed Brown 
2376c924f48SJed Brown #undef __FUNCT__
23869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
23969a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2400971522cSBarry Smith {
2410971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2420971522cSBarry Smith   PetscErrorCode    ierr;
2435a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2447287d2a3SDmitry Karpeev   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE;
2456c924f48SJed Brown   PetscInt          i;
2460971522cSBarry Smith 
2470971522cSBarry Smith   PetscFunctionBegin;
2487287d2a3SDmitry Karpeev   /*
2497287d2a3SDmitry Karpeev    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
2507287d2a3SDmitry Karpeev    Should probably be rewritten.
2517287d2a3SDmitry Karpeev    */
2527287d2a3SDmitry Karpeev   if (!ilink || jac->reset) {
253d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
254d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
255bafc1b83SMatthew G Knepley       PetscInt     numFields, f, i, j;
2560784a22cSJed Brown       char       **fieldNames;
2577b62db95SJungho Lee       IS          *fields;
258e7c4fc90SDmitry Karpeev       DM          *dms;
259bafc1b83SMatthew G Knepley       DM           subdm[128];
260bafc1b83SMatthew G Knepley       PetscBool    flg;
261bafc1b83SMatthew G Knepley 
26216621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
263bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
264bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE; ; i++) {
265bafc1b83SMatthew G Knepley         PetscInt ifields[128];
266bafc1b83SMatthew G Knepley         IS       compField;
267bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
268bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
269bafc1b83SMatthew G Knepley 
270*8caf3d72SBarry Smith         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
271bafc1b83SMatthew G Knepley         ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
272bafc1b83SMatthew G Knepley         if (!flg) break;
273bafc1b83SMatthew G Knepley         if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
274bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
275bafc1b83SMatthew G Knepley         if (nfields == 1) {
276bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
2774d2389d7SMatthew G Knepley           /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr);
2784d2389d7SMatthew G Knepley              ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */
279bafc1b83SMatthew G Knepley         } else {
280*8caf3d72SBarry Smith           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
281bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
2824d2389d7SMatthew G Knepley           /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr);
2834d2389d7SMatthew G Knepley              ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */
2847287d2a3SDmitry Karpeev         }
285bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
286bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
287bafc1b83SMatthew G Knepley           f = ifields[j];
2887b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
2897b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
2907b62db95SJungho Lee         }
291bafc1b83SMatthew G Knepley       }
292bafc1b83SMatthew G Knepley       if (i == 0) {
293bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
294bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
295bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
296bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
297bafc1b83SMatthew G Knepley         }
298bafc1b83SMatthew G Knepley       } else {
299bafc1b83SMatthew G Knepley         ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr);
300bafc1b83SMatthew G Knepley         for (j = 0; j < i; ++j) {
301bafc1b83SMatthew G Knepley           dms[j] = subdm[j];
302bafc1b83SMatthew G Knepley         }
303bafc1b83SMatthew G Knepley       }
3047b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
3057b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
306e7c4fc90SDmitry Karpeev       if (dms) {
3078b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
308bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
3097287d2a3SDmitry Karpeev           const char *prefix;
3107287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr);
3117287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);     CHKERRQ(ierr);
3127b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
3137b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
3147287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr);
315e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]); CHKERRQ(ierr);
3162fa5ba8aSJed Brown         }
3177b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
3188b8307b2SJed Brown       }
31966ffff09SJed Brown     } else {
320521d7252SBarry Smith       if (jac->bs <= 0) {
321704ba839SBarry Smith         if (pc->pmat) {
322521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
323704ba839SBarry Smith         } else {
324704ba839SBarry Smith           jac->bs = 1;
325704ba839SBarry Smith         }
326521d7252SBarry Smith       }
327d32f9abdSBarry Smith 
3287287d2a3SDmitry Karpeev       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr);
3296ce1633cSBarry Smith       if (stokes) {
3306ce1633cSBarry Smith         IS       zerodiags,rest;
3316ce1633cSBarry Smith         PetscInt nmin,nmax;
3326ce1633cSBarry Smith 
3336ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
3346ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
3356ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
3367287d2a3SDmitry Karpeev         if (jac->reset) {
3377287d2a3SDmitry Karpeev           jac->head->is       = rest;
3387287d2a3SDmitry Karpeev           jac->head->next->is = zerodiags;
3397287d2a3SDmitry Karpeev         }
3407287d2a3SDmitry Karpeev         else {
3416ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
3426ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
3437287d2a3SDmitry Karpeev         }
344fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
345fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
3466ce1633cSBarry Smith       } else {
3477287d2a3SDmitry Karpeev         if (jac->reset)
3487287d2a3SDmitry Karpeev           SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
3497287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
350d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
351d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
3526c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
3536c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
354d32f9abdSBarry Smith         }
3557287d2a3SDmitry Karpeev         if (fieldsplit_default || !jac->splitdefined) {
356d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
357db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
3586c924f48SJed Brown             char splitname[8];
359*8caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
3605d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
36179416396SBarry Smith           }
3625d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
363521d7252SBarry Smith         }
36466ffff09SJed Brown       }
3656ce1633cSBarry Smith     }
366edf189efSBarry Smith   } else if (jac->nsplits == 1) {
367edf189efSBarry Smith     if (ilink->is) {
368edf189efSBarry Smith       IS       is2;
369edf189efSBarry Smith       PetscInt nmin,nmax;
370edf189efSBarry Smith 
371edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
372edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
373db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
374fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
3757b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
376edf189efSBarry Smith   }
377d0af7cd3SBarry Smith 
378d0af7cd3SBarry Smith 
3797b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
38069a612a9SBarry Smith   PetscFunctionReturn(0);
38169a612a9SBarry Smith }
38269a612a9SBarry Smith 
383514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg);
384514bf10dSMatthew G Knepley 
38569a612a9SBarry Smith #undef __FUNCT__
38669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
38769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
38869a612a9SBarry Smith {
38969a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
39069a612a9SBarry Smith   PetscErrorCode    ierr;
3915a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
3922c9966d7SBarry Smith   PetscInt          i,nsplit;
39369a612a9SBarry Smith   MatStructure      flag = pc->flag;
3945f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
39569a612a9SBarry Smith 
39669a612a9SBarry Smith   PetscFunctionBegin;
39769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
39897bbdb24SBarry Smith   nsplit = jac->nsplits;
3995a9f2f41SSatish Balay   ilink  = jac->head;
40097bbdb24SBarry Smith 
40197bbdb24SBarry Smith   /* get the matrices for each split */
402704ba839SBarry Smith   if (!jac->issetup) {
4031b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
40497bbdb24SBarry Smith 
405704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
406704ba839SBarry Smith 
4075d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4082c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4092c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
4102c9966d7SBarry Smith     }
41151f519a2SBarry Smith     bs     = jac->bs;
41297bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
4131b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4141b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4151b9fc7fcSBarry Smith       if (jac->defaultsplit) {
416704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
4175f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
418704ba839SBarry Smith       } else if (!ilink->is) {
419ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4205f4ab4e1SJungho Lee           PetscInt   *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
4215a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
4225f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
4231b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
4241b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
4251b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
4265f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
42797bbdb24SBarry Smith             }
42897bbdb24SBarry Smith           }
429d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
4305f4ab4e1SJungho Lee           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
431ccb205f8SBarry Smith         } else {
432704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
4335f4ab4e1SJungho Lee           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
434ccb205f8SBarry Smith        }
4353e197d65SBarry Smith       }
436704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
4375f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
4385f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
439704ba839SBarry Smith       ilink = ilink->next;
4401b9fc7fcSBarry Smith     }
4411b9fc7fcSBarry Smith   }
4421b9fc7fcSBarry Smith 
443704ba839SBarry Smith   ilink  = jac->head;
44497bbdb24SBarry Smith   if (!jac->pmat) {
445bdddcaaaSMatthew G Knepley     Vec xtmp;
446bdddcaaaSMatthew G Knepley 
447bdddcaaaSMatthew G Knepley     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
448cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
449bdddcaaaSMatthew G Knepley     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
450cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
451bdddcaaaSMatthew G Knepley       MatNullSpace sp;
452bdddcaaaSMatthew G Knepley 
453a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
454a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr);
455a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
456a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
457a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
458a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
459a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
460a3df900dSMatthew G Knepley         }
461a3df900dSMatthew G Knepley       } else {
4625f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
463a3df900dSMatthew G Knepley       }
464bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
465bdddcaaaSMatthew G Knepley       ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
466bdddcaaaSMatthew G Knepley       ilink->x = jac->x[i]; ilink->y = jac->y[i];
467bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
468bdddcaaaSMatthew G Knepley       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
469ed1f0337SMatthew G Knepley       /* Check for null space attached to IS */
470bafc1b83SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr);
471bafc1b83SMatthew G Knepley       if (sp) {
472bafc1b83SMatthew G Knepley         ierr  = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
473bafc1b83SMatthew G Knepley       }
474ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr);
475ed1f0337SMatthew G Knepley       if (sp) {
476ed1f0337SMatthew G Knepley         ierr  = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
477ed1f0337SMatthew G Knepley       }
478704ba839SBarry Smith       ilink = ilink->next;
479cf502942SBarry Smith     }
480bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
48197bbdb24SBarry Smith   } else {
482cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
483a3df900dSMatthew G Knepley       Mat pmat;
484a3df900dSMatthew G Knepley 
485a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
486a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr);
487a3df900dSMatthew G Knepley       if (!pmat) {
4885f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
489a3df900dSMatthew G Knepley       }
490704ba839SBarry Smith       ilink = ilink->next;
491cf502942SBarry Smith     }
49297bbdb24SBarry Smith   }
493519d70e2SJed Brown   if (jac->realdiagonal) {
494519d70e2SJed Brown     ilink = jac->head;
495519d70e2SJed Brown     if (!jac->mat) {
496519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
497519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4985f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
499519d70e2SJed Brown         ilink = ilink->next;
500519d70e2SJed Brown       }
501519d70e2SJed Brown     } else {
502519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5035f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
504519d70e2SJed Brown         ilink = ilink->next;
505519d70e2SJed Brown       }
506519d70e2SJed Brown     }
507519d70e2SJed Brown   } else {
508519d70e2SJed Brown     jac->mat = jac->pmat;
509519d70e2SJed Brown   }
51097bbdb24SBarry Smith 
5116c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
51268dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
51368dd23aaSBarry Smith     ilink  = jac->head;
51468dd23aaSBarry Smith     if (!jac->Afield) {
51568dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
51668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5174aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
51868dd23aaSBarry Smith         ilink = ilink->next;
51968dd23aaSBarry Smith       }
52068dd23aaSBarry Smith     } else {
52168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5224aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
52368dd23aaSBarry Smith         ilink = ilink->next;
52468dd23aaSBarry Smith       }
52568dd23aaSBarry Smith     }
52668dd23aaSBarry Smith   }
52768dd23aaSBarry Smith 
5283b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
5293b224e63SBarry Smith     IS       ccis;
5304aa3045dSJed Brown     PetscInt rstart,rend;
531093c86ffSJed Brown     char     lscname[256];
532093c86ffSJed Brown     PetscObject LSC_L;
533e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
53468dd23aaSBarry Smith 
535e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
536e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
537e6cab6aaSJed Brown 
5383b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
5393b224e63SBarry Smith     if (jac->schur) {
5403b224e63SBarry Smith       ilink = jac->head;
54149bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5424aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
543fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5443b224e63SBarry Smith       ilink = ilink->next;
54549bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5464aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
547fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
548a3df900dSMatthew G Knepley       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr);
549084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
5503b224e63SBarry Smith      } else {
5511cee3971SBarry Smith       KSP ksp;
552bafc1b83SMatthew G Knepley       const char  *Dprefix;
553bafc1b83SMatthew G Knepley       char schurprefix[256];
554514bf10dSMatthew G Knepley       char schurtestoption[256];
555bdddcaaaSMatthew G Knepley       MatNullSpace sp;
556514bf10dSMatthew G Knepley       PetscBool    flg;
5573b224e63SBarry Smith 
558a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
5593b224e63SBarry Smith       ilink = jac->head;
56049bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5614aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
562fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5633b224e63SBarry Smith       ilink = ilink->next;
56449bb4cd7SJungho Lee       ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5654aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
566fcfd50ebSBarry Smith       ierr = ISDestroy(&ccis);CHKERRQ(ierr);
56720252d06SBarry Smith 
56820252d06SBarry Smith       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */
56920252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
57020252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
57120252d06SBarry Smith       ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr);
57220252d06SBarry Smith       ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
57320252d06SBarry Smith       /* Indent this deeper to emphasize the "inner" nature of this solver. */
57420252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr);
57520252d06SBarry Smith       ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr);
57620252d06SBarry Smith       ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
57720252d06SBarry Smith 
578bdddcaaaSMatthew G Knepley       ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr);
57920252d06SBarry Smith       if (sp) {
58020252d06SBarry Smith         ierr  = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
58120252d06SBarry Smith       }
58220252d06SBarry Smith 
58320252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
584514bf10dSMatthew G Knepley       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr);
585514bf10dSMatthew G Knepley       if (flg) {
586514bf10dSMatthew G Knepley         DM dmInner;
587514bf10dSMatthew G Knepley 
588514bf10dSMatthew G Knepley         /* Set DM for new solver */
589bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
590bafc1b83SMatthew G Knepley         ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr);
5917287d2a3SDmitry Karpeev         ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
592514bf10dSMatthew G Knepley       } else {
593514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr);
594514bf10dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
595bafc1b83SMatthew G Knepley       }
59620b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
5973b224e63SBarry Smith 
5983b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);               CHKERRQ(ierr);
5999005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
60020252d06SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);           CHKERRQ(ierr);
601084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
602084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
6037233a360SDmitry Karpeev         PC pcschur;
6047233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
6057233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
606084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
607e69d4d44SBarry Smith       }
6087287d2a3SDmitry Karpeev       ierr  = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr);
6097287d2a3SDmitry Karpeev       ierr  = KSPSetOptionsPrefix(jac->kspschur,         Dprefix); CHKERRQ(ierr);
6103b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
61120b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
61220b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
6133b224e63SBarry Smith     }
614093c86ffSJed Brown 
615093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
616*8caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
617093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
618093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
619093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
620*8caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
621093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
622093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
623093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
6243b224e63SBarry Smith   } else {
62568bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
62697bbdb24SBarry Smith     i    = 0;
6275a9f2f41SSatish Balay     ilink = jac->head;
6285a9f2f41SSatish Balay     while (ilink) {
629519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
6303b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
631c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
63297bbdb24SBarry Smith       i++;
6335a9f2f41SSatish Balay       ilink = ilink->next;
6340971522cSBarry Smith     }
6353b224e63SBarry Smith   }
6363b224e63SBarry Smith 
637c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
6380971522cSBarry Smith   PetscFunctionReturn(0);
6390971522cSBarry Smith }
6400971522cSBarry Smith 
6415a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
642ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
643ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
6445a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
645ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
646ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
64779416396SBarry Smith 
6480971522cSBarry Smith #undef __FUNCT__
6493b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
6503b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
6513b224e63SBarry Smith {
6523b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6533b224e63SBarry Smith   PetscErrorCode    ierr;
6543b224e63SBarry Smith   KSP               ksp;
6553b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
6563b224e63SBarry Smith 
6573b224e63SBarry Smith   PetscFunctionBegin;
6583b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
6593b224e63SBarry Smith 
660c5d2311dSJed Brown   switch (jac->schurfactorization) {
661c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
662a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
663c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
664c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
665c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
666c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
667c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
668c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
669c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
670c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
671c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
672c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
673c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
674c5d2311dSJed Brown       break;
675c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
676a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
677c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
678c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
679c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
680c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
681c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
682c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
683c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
684c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
685c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
686c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
687c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
688c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
689c5d2311dSJed Brown       break;
690c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
691a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
692c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
693c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
694c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
695c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
696c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
697c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
698c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
699c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
700c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
701c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
702c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
703c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
704c5d2311dSJed Brown       break;
705c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_FULL:
706a04f6461SBarry 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 */
7073b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7083b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7093b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
7103b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
7113b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
7123b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7133b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7143b224e63SBarry Smith 
7153b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
7163b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7173b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7183b224e63SBarry Smith 
7193b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
7203b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
7213b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
7223b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7233b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
724c5d2311dSJed Brown   }
7253b224e63SBarry Smith   PetscFunctionReturn(0);
7263b224e63SBarry Smith }
7273b224e63SBarry Smith 
7283b224e63SBarry Smith #undef __FUNCT__
7290971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
7300971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
7310971522cSBarry Smith {
7320971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7330971522cSBarry Smith   PetscErrorCode    ierr;
7345a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
735939b8a20SBarry Smith   PetscInt          cnt,bs;
7360971522cSBarry Smith 
7370971522cSBarry Smith   PetscFunctionBegin;
7384442daceSBarry Smith   x->map->bs = jac->bs;
7394442daceSBarry Smith   y->map->bs = jac->bs;
74051f519a2SBarry Smith   CHKMEMQ;
74179416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
7421b9fc7fcSBarry Smith     if (jac->defaultsplit) {
743939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
744939b8a20SBarry 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);
745939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
746939b8a20SBarry 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);
7470971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
7485a9f2f41SSatish Balay       while (ilink) {
7495a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7505a9f2f41SSatish Balay         ilink = ilink->next;
7510971522cSBarry Smith       }
7520971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
7531b9fc7fcSBarry Smith     } else {
754efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
7555a9f2f41SSatish Balay       while (ilink) {
7565a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7575a9f2f41SSatish Balay         ilink = ilink->next;
7581b9fc7fcSBarry Smith       }
7591b9fc7fcSBarry Smith     }
76016913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
76179416396SBarry Smith     if (!jac->w1) {
76279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
76379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
76479416396SBarry Smith     }
765efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
7665a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7673e197d65SBarry Smith     cnt = 1;
7685a9f2f41SSatish Balay     while (ilink->next) {
7695a9f2f41SSatish Balay       ilink = ilink->next;
7703e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7713e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7723e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7733e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7743e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7753e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7763e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7773e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7783e197d65SBarry Smith     }
77951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
78011755939SBarry Smith       cnt -= 2;
78151f519a2SBarry Smith       while (ilink->previous) {
78251f519a2SBarry Smith         ilink = ilink->previous;
78311755939SBarry Smith         /* compute the residual only over the part of the vector needed */
78411755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
78511755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
78611755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
78711755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
78811755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
78911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
79011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
79151f519a2SBarry Smith       }
79211755939SBarry Smith     }
79365e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
79451f519a2SBarry Smith   CHKMEMQ;
7950971522cSBarry Smith   PetscFunctionReturn(0);
7960971522cSBarry Smith }
7970971522cSBarry Smith 
798421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
799ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
800ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
801421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
802ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
803ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
804421e10b8SBarry Smith 
805421e10b8SBarry Smith #undef __FUNCT__
8068c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
807421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
808421e10b8SBarry Smith {
809421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
810421e10b8SBarry Smith   PetscErrorCode    ierr;
811421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
812939b8a20SBarry Smith   PetscInt          bs;
813421e10b8SBarry Smith 
814421e10b8SBarry Smith   PetscFunctionBegin;
815421e10b8SBarry Smith   CHKMEMQ;
816421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
817421e10b8SBarry Smith     if (jac->defaultsplit) {
818939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
819939b8a20SBarry 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);
820939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
821939b8a20SBarry 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);
822421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
823421e10b8SBarry Smith       while (ilink) {
824421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
825421e10b8SBarry Smith 	ilink = ilink->next;
826421e10b8SBarry Smith       }
827421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
828421e10b8SBarry Smith     } else {
829421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
830421e10b8SBarry Smith       while (ilink) {
831421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
832421e10b8SBarry Smith 	ilink = ilink->next;
833421e10b8SBarry Smith       }
834421e10b8SBarry Smith     }
835421e10b8SBarry Smith   } else {
836421e10b8SBarry Smith     if (!jac->w1) {
837421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
838421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
839421e10b8SBarry Smith     }
840421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
841421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
842421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
843421e10b8SBarry Smith       while (ilink->next) {
844421e10b8SBarry Smith         ilink = ilink->next;
8459989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
846421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
847421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
848421e10b8SBarry Smith       }
849421e10b8SBarry Smith       while (ilink->previous) {
850421e10b8SBarry Smith         ilink = ilink->previous;
8519989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
852421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
853421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
854421e10b8SBarry Smith       }
855421e10b8SBarry Smith     } else {
856421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
857421e10b8SBarry Smith 	ilink = ilink->next;
858421e10b8SBarry Smith       }
859421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
860421e10b8SBarry Smith       while (ilink->previous) {
861421e10b8SBarry Smith 	ilink = ilink->previous;
8629989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
863421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
864421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
865421e10b8SBarry Smith       }
866421e10b8SBarry Smith     }
867421e10b8SBarry Smith   }
868421e10b8SBarry Smith   CHKMEMQ;
869421e10b8SBarry Smith   PetscFunctionReturn(0);
870421e10b8SBarry Smith }
871421e10b8SBarry Smith 
8720971522cSBarry Smith #undef __FUNCT__
873574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
874574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8750971522cSBarry Smith {
8760971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8770971522cSBarry Smith   PetscErrorCode    ierr;
8785a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8790971522cSBarry Smith 
8800971522cSBarry Smith   PetscFunctionBegin;
8815a9f2f41SSatish Balay   while (ilink) {
882574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
883fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
884fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
885fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
886fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
887b5787286SJed Brown     ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
8885a9f2f41SSatish Balay     next = ilink->next;
8895a9f2f41SSatish Balay     ilink = next;
8900971522cSBarry Smith   }
89105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
892574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
893574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
894574deadeSBarry Smith   } else if (jac->mat) {
895574deadeSBarry Smith     jac->mat = PETSC_NULL;
896574deadeSBarry Smith   }
89797bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
89868dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8996bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
9006bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
9016bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
9026bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
9036bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
9046bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
9056bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
90663ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
907574deadeSBarry Smith   PetscFunctionReturn(0);
908574deadeSBarry Smith }
909574deadeSBarry Smith 
910574deadeSBarry Smith #undef __FUNCT__
911574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
912574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
913574deadeSBarry Smith {
914574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
915574deadeSBarry Smith   PetscErrorCode    ierr;
916574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
917574deadeSBarry Smith 
918574deadeSBarry Smith   PetscFunctionBegin;
919574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
920574deadeSBarry Smith   while (ilink) {
9216bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
922574deadeSBarry Smith     next = ilink->next;
923574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
924574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
9255d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
926574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
927574deadeSBarry Smith     ilink = next;
928574deadeSBarry Smith   }
929574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
930c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
9317b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
9327b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
9337b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
9347b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
9357b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
936ab1df9f5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr);
937c9c6ffaaSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr);
9380971522cSBarry Smith   PetscFunctionReturn(0);
9390971522cSBarry Smith }
9400971522cSBarry Smith 
9410971522cSBarry Smith #undef __FUNCT__
9420971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
9430971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
9440971522cSBarry Smith {
9451b9fc7fcSBarry Smith   PetscErrorCode  ierr;
9466c924f48SJed Brown   PetscInt        bs;
947bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
9489dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
9493b224e63SBarry Smith   PCCompositeType ctype;
950e7c4fc90SDmitry Karpeev   DM              ddm;
951e7c4fc90SDmitry Karpeev   char            ddm_name[1025];
9521b9fc7fcSBarry Smith 
9530971522cSBarry Smith   PetscFunctionBegin;
9541b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
955e7c4fc90SDmitry Karpeev   if (pc->dm) {
956e7c4fc90SDmitry Karpeev     /* Allow the user to request a decomposition DM by name */
957e7c4fc90SDmitry Karpeev     ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
958731c8d9eSDmitry Karpeev     ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
959e7c4fc90SDmitry Karpeev     if (flg) {
96016621825SDmitry Karpeev       ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
961e7c4fc90SDmitry Karpeev       if (!ddm) {
96216621825SDmitry Karpeev         SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name);
963e7c4fc90SDmitry Karpeev       }
96416621825SDmitry Karpeev       ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr);
965e7c4fc90SDmitry Karpeev       ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
966e7c4fc90SDmitry Karpeev     }
967e7c4fc90SDmitry Karpeev   }
968acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
96951f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
97051f519a2SBarry Smith   if (flg) {
97151f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
97251f519a2SBarry Smith   }
973704ba839SBarry Smith 
974435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
975c0adfefeSBarry Smith   if (stokes) {
976c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
977c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
978c0adfefeSBarry Smith   }
979c0adfefeSBarry Smith 
9803b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
9813b224e63SBarry Smith   if (flg) {
9823b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
9833b224e63SBarry Smith   }
984c30613efSMatthew Knepley   /* Only setup fields once */
985c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
986d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
987d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
9886c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9896c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
990d32f9abdSBarry Smith   }
991c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
992c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
993c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
994c9c6ffaaSJed 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);
995084e4875SJed 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);
996c5d2311dSJed Brown   }
9971b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9980971522cSBarry Smith   PetscFunctionReturn(0);
9990971522cSBarry Smith }
10000971522cSBarry Smith 
10010971522cSBarry Smith /*------------------------------------------------------------------------------------*/
10020971522cSBarry Smith 
10030971522cSBarry Smith EXTERN_C_BEGIN
10040971522cSBarry Smith #undef __FUNCT__
10050971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
10065d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
10070971522cSBarry Smith {
100897bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10090971522cSBarry Smith   PetscErrorCode    ierr;
10105a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
101169a612a9SBarry Smith   char              prefix[128];
10125d4c12cdSJungho Lee   PetscInt          i;
10130971522cSBarry Smith 
10140971522cSBarry Smith   PetscFunctionBegin;
10156c924f48SJed Brown   if (jac->splitdefined) {
10166c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10176c924f48SJed Brown     PetscFunctionReturn(0);
10186c924f48SJed Brown   }
101951f519a2SBarry Smith   for (i=0; i<n; i++) {
1020e32f2f54SBarry 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);
1021e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
102251f519a2SBarry Smith   }
1023704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1024a04f6461SBarry Smith   if (splitname) {
1025db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1026a04f6461SBarry Smith   } else {
1027a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1028a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1029a04f6461SBarry Smith   }
1030704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
10315a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
10325d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
10335d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
10345a9f2f41SSatish Balay   ilink->nfields = n;
10355a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
10367adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
103720252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
10385a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10399005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
104069a612a9SBarry Smith 
1041*8caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
10425a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
10430971522cSBarry Smith 
10440971522cSBarry Smith   if (!next) {
10455a9f2f41SSatish Balay     jac->head       = ilink;
104651f519a2SBarry Smith     ilink->previous = PETSC_NULL;
10470971522cSBarry Smith   } else {
10480971522cSBarry Smith     while (next->next) {
10490971522cSBarry Smith       next = next->next;
10500971522cSBarry Smith     }
10515a9f2f41SSatish Balay     next->next      = ilink;
105251f519a2SBarry Smith     ilink->previous = next;
10530971522cSBarry Smith   }
10540971522cSBarry Smith   jac->nsplits++;
10550971522cSBarry Smith   PetscFunctionReturn(0);
10560971522cSBarry Smith }
10570971522cSBarry Smith EXTERN_C_END
10580971522cSBarry Smith 
1059e69d4d44SBarry Smith EXTERN_C_BEGIN
1060e69d4d44SBarry Smith #undef __FUNCT__
1061e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
10627087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1063e69d4d44SBarry Smith {
1064e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1065e69d4d44SBarry Smith   PetscErrorCode ierr;
1066e69d4d44SBarry Smith 
1067e69d4d44SBarry Smith   PetscFunctionBegin;
1068519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1069e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1070e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
107113e0d083SBarry Smith   if (n) *n = jac->nsplits;
1072e69d4d44SBarry Smith   PetscFunctionReturn(0);
1073e69d4d44SBarry Smith }
1074e69d4d44SBarry Smith EXTERN_C_END
10750971522cSBarry Smith 
10760971522cSBarry Smith EXTERN_C_BEGIN
10770971522cSBarry Smith #undef __FUNCT__
107869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
10797087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
10800971522cSBarry Smith {
10810971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10820971522cSBarry Smith   PetscErrorCode    ierr;
10830971522cSBarry Smith   PetscInt          cnt = 0;
10845a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
10850971522cSBarry Smith 
10860971522cSBarry Smith   PetscFunctionBegin;
10875d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
10885a9f2f41SSatish Balay   while (ilink) {
10895a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
10905a9f2f41SSatish Balay     ilink = ilink->next;
10910971522cSBarry Smith   }
10925d480477SMatthew 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);
109313e0d083SBarry Smith   if (n) *n = jac->nsplits;
10940971522cSBarry Smith   PetscFunctionReturn(0);
10950971522cSBarry Smith }
10960971522cSBarry Smith EXTERN_C_END
10970971522cSBarry Smith 
1098704ba839SBarry Smith EXTERN_C_BEGIN
1099704ba839SBarry Smith #undef __FUNCT__
1100704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
11017087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1102704ba839SBarry Smith {
1103704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1104704ba839SBarry Smith   PetscErrorCode    ierr;
1105704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1106704ba839SBarry Smith   char              prefix[128];
1107704ba839SBarry Smith 
1108704ba839SBarry Smith   PetscFunctionBegin;
11096c924f48SJed Brown   if (jac->splitdefined) {
11106c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
11116c924f48SJed Brown     PetscFunctionReturn(0);
11126c924f48SJed Brown   }
111316913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1114a04f6461SBarry Smith   if (splitname) {
1115db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1116a04f6461SBarry Smith   } else {
1117b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1118b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1119a04f6461SBarry Smith   }
11201ab39975SBarry Smith   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1121b5787286SJed Brown   ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
1122b5787286SJed Brown   ilink->is      = is;
1123b5787286SJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1124b5787286SJed Brown   ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1125b5787286SJed Brown   ilink->is_col  = is;
1126704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1127704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
112820252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1129704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
11309005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1131704ba839SBarry Smith 
1132*8caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1133704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1134704ba839SBarry Smith 
1135704ba839SBarry Smith   if (!next) {
1136704ba839SBarry Smith     jac->head       = ilink;
1137704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1138704ba839SBarry Smith   } else {
1139704ba839SBarry Smith     while (next->next) {
1140704ba839SBarry Smith       next = next->next;
1141704ba839SBarry Smith     }
1142704ba839SBarry Smith     next->next      = ilink;
1143704ba839SBarry Smith     ilink->previous = next;
1144704ba839SBarry Smith   }
1145704ba839SBarry Smith   jac->nsplits++;
1146704ba839SBarry Smith 
1147704ba839SBarry Smith   PetscFunctionReturn(0);
1148704ba839SBarry Smith }
1149704ba839SBarry Smith EXTERN_C_END
1150704ba839SBarry Smith 
11510971522cSBarry Smith #undef __FUNCT__
11520971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
11530971522cSBarry Smith /*@
11540971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
11550971522cSBarry Smith 
1156ad4df100SBarry Smith     Logically Collective on PC
11570971522cSBarry Smith 
11580971522cSBarry Smith     Input Parameters:
11590971522cSBarry Smith +   pc  - the preconditioner context
1160a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
11610971522cSBarry Smith .   n - the number of fields in this split
1162db4c96c1SJed Brown -   fields - the fields in this split
11630971522cSBarry Smith 
11640971522cSBarry Smith     Level: intermediate
11650971522cSBarry Smith 
1166d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1167d32f9abdSBarry Smith 
11687287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1169d32f9abdSBarry 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
1170d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1171d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1172d32f9abdSBarry Smith 
1173db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1174db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1175db4c96c1SJed Brown 
11765d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
11775d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
11785d4c12cdSJungho Lee      available when this routine is called.
11795d4c12cdSJungho Lee 
1180d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
11810971522cSBarry Smith 
11820971522cSBarry Smith @*/
11835d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
11840971522cSBarry Smith {
11854ac538c5SBarry Smith   PetscErrorCode ierr;
11860971522cSBarry Smith 
11870971522cSBarry Smith   PetscFunctionBegin;
11880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1189db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1190db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1191db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
11925d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
11930971522cSBarry Smith   PetscFunctionReturn(0);
11940971522cSBarry Smith }
11950971522cSBarry Smith 
11960971522cSBarry Smith #undef __FUNCT__
1197704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1198704ba839SBarry Smith /*@
1199704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1200704ba839SBarry Smith 
1201ad4df100SBarry Smith     Logically Collective on PC
1202704ba839SBarry Smith 
1203704ba839SBarry Smith     Input Parameters:
1204704ba839SBarry Smith +   pc  - the preconditioner context
1205a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1206db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1207704ba839SBarry Smith 
1208d32f9abdSBarry Smith 
1209a6ffb8dbSJed Brown     Notes:
1210a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1211a6ffb8dbSJed Brown 
1212db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1213db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1214d32f9abdSBarry Smith 
1215704ba839SBarry Smith     Level: intermediate
1216704ba839SBarry Smith 
1217704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1218704ba839SBarry Smith 
1219704ba839SBarry Smith @*/
12207087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1221704ba839SBarry Smith {
12224ac538c5SBarry Smith   PetscErrorCode ierr;
1223704ba839SBarry Smith 
1224704ba839SBarry Smith   PetscFunctionBegin;
12250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12267b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1227db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
12284ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1229704ba839SBarry Smith   PetscFunctionReturn(0);
1230704ba839SBarry Smith }
1231704ba839SBarry Smith 
1232704ba839SBarry Smith #undef __FUNCT__
123357a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
123457a9adfeSMatthew G Knepley /*@
123557a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
123657a9adfeSMatthew G Knepley 
123757a9adfeSMatthew G Knepley     Logically Collective on PC
123857a9adfeSMatthew G Knepley 
123957a9adfeSMatthew G Knepley     Input Parameters:
124057a9adfeSMatthew G Knepley +   pc  - the preconditioner context
124157a9adfeSMatthew G Knepley -   splitname - name of this split
124257a9adfeSMatthew G Knepley 
124357a9adfeSMatthew G Knepley     Output Parameter:
124457a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
124557a9adfeSMatthew G Knepley 
124657a9adfeSMatthew G Knepley     Level: intermediate
124757a9adfeSMatthew G Knepley 
124857a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
124957a9adfeSMatthew G Knepley 
125057a9adfeSMatthew G Knepley @*/
125157a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
125257a9adfeSMatthew G Knepley {
125357a9adfeSMatthew G Knepley   PetscErrorCode ierr;
125457a9adfeSMatthew G Knepley 
125557a9adfeSMatthew G Knepley   PetscFunctionBegin;
125657a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
125757a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
125857a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
125957a9adfeSMatthew G Knepley   {
126057a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
126157a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
126257a9adfeSMatthew G Knepley     PetscBool         found;
126357a9adfeSMatthew G Knepley 
126457a9adfeSMatthew G Knepley     *is = PETSC_NULL;
126557a9adfeSMatthew G Knepley     while(ilink) {
126657a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
126757a9adfeSMatthew G Knepley       if (found) {
126857a9adfeSMatthew G Knepley         *is = ilink->is;
126957a9adfeSMatthew G Knepley         break;
127057a9adfeSMatthew G Knepley       }
127157a9adfeSMatthew G Knepley       ilink = ilink->next;
127257a9adfeSMatthew G Knepley     }
127357a9adfeSMatthew G Knepley   }
127457a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
127557a9adfeSMatthew G Knepley }
127657a9adfeSMatthew G Knepley 
127757a9adfeSMatthew G Knepley #undef __FUNCT__
127851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
127951f519a2SBarry Smith /*@
128051f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
128151f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
128251f519a2SBarry Smith 
1283ad4df100SBarry Smith     Logically Collective on PC
128451f519a2SBarry Smith 
128551f519a2SBarry Smith     Input Parameters:
128651f519a2SBarry Smith +   pc  - the preconditioner context
128751f519a2SBarry Smith -   bs - the block size
128851f519a2SBarry Smith 
128951f519a2SBarry Smith     Level: intermediate
129051f519a2SBarry Smith 
129151f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
129251f519a2SBarry Smith 
129351f519a2SBarry Smith @*/
12947087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
129551f519a2SBarry Smith {
12964ac538c5SBarry Smith   PetscErrorCode ierr;
129751f519a2SBarry Smith 
129851f519a2SBarry Smith   PetscFunctionBegin;
12990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1300c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
13014ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
130251f519a2SBarry Smith   PetscFunctionReturn(0);
130351f519a2SBarry Smith }
130451f519a2SBarry Smith 
130551f519a2SBarry Smith #undef __FUNCT__
130669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
13070971522cSBarry Smith /*@C
130869a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
13090971522cSBarry Smith 
131069a612a9SBarry Smith    Collective on KSP
13110971522cSBarry Smith 
13120971522cSBarry Smith    Input Parameter:
13130971522cSBarry Smith .  pc - the preconditioner context
13140971522cSBarry Smith 
13150971522cSBarry Smith    Output Parameters:
131613e0d083SBarry Smith +  n - the number of splits
131769a612a9SBarry Smith -  pc - the array of KSP contexts
13180971522cSBarry Smith 
13190971522cSBarry Smith    Note:
1320d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1321d32f9abdSBarry Smith    (not the KSP just the array that contains them).
13220971522cSBarry Smith 
132369a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
13240971522cSBarry Smith 
1325196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1326196cc216SBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the
1327196cc216SBarry Smith       KSP array must be.
1328196cc216SBarry Smith 
1329196cc216SBarry Smith 
13300971522cSBarry Smith    Level: advanced
13310971522cSBarry Smith 
13320971522cSBarry Smith .seealso: PCFIELDSPLIT
13330971522cSBarry Smith @*/
13347087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
13350971522cSBarry Smith {
13364ac538c5SBarry Smith   PetscErrorCode ierr;
13370971522cSBarry Smith 
13380971522cSBarry Smith   PetscFunctionBegin;
13390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
134013e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
13414ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
13420971522cSBarry Smith   PetscFunctionReturn(0);
13430971522cSBarry Smith }
13440971522cSBarry Smith 
1345e69d4d44SBarry Smith #undef __FUNCT__
1346e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1347e69d4d44SBarry Smith /*@
1348e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1349a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1350e69d4d44SBarry Smith 
1351e69d4d44SBarry Smith     Collective on PC
1352e69d4d44SBarry Smith 
1353e69d4d44SBarry Smith     Input Parameters:
1354e69d4d44SBarry Smith +   pc  - the preconditioner context
1355fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1356084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1357084e4875SJed Brown 
1358e69d4d44SBarry Smith     Options Database:
1359084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1360e69d4d44SBarry Smith 
1361fd1303e9SJungho Lee     Notes:
1362fd1303e9SJungho Lee $    If ptype is
1363fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1364fd1303e9SJungho Lee $             to this function).
1365fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1366fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1367fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1368fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1369fd1303e9SJungho Lee $             preconditioner
1370fd1303e9SJungho Lee 
1371fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1372fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1373fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1374fd1303e9SJungho Lee 
1375fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1376fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1377fd1303e9SJungho Lee 
1378e69d4d44SBarry Smith     Level: intermediate
1379e69d4d44SBarry Smith 
1380fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1381e69d4d44SBarry Smith 
1382e69d4d44SBarry Smith @*/
13837087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1384e69d4d44SBarry Smith {
13854ac538c5SBarry Smith   PetscErrorCode ierr;
1386e69d4d44SBarry Smith 
1387e69d4d44SBarry Smith   PetscFunctionBegin;
13880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13894ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1390e69d4d44SBarry Smith   PetscFunctionReturn(0);
1391e69d4d44SBarry Smith }
1392e69d4d44SBarry Smith 
1393e69d4d44SBarry Smith EXTERN_C_BEGIN
1394e69d4d44SBarry Smith #undef __FUNCT__
1395e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
13967087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1397e69d4d44SBarry Smith {
1398e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1399084e4875SJed Brown   PetscErrorCode  ierr;
1400e69d4d44SBarry Smith 
1401e69d4d44SBarry Smith   PetscFunctionBegin;
1402084e4875SJed Brown   jac->schurpre = ptype;
1403084e4875SJed Brown   if (pre) {
14046bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1405084e4875SJed Brown     jac->schur_user = pre;
1406084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1407084e4875SJed Brown   }
1408e69d4d44SBarry Smith   PetscFunctionReturn(0);
1409e69d4d44SBarry Smith }
1410e69d4d44SBarry Smith EXTERN_C_END
1411e69d4d44SBarry Smith 
141230ad9308SMatthew Knepley #undef __FUNCT__
1413c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1414ab1df9f5SJed Brown /*@
1415c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1416ab1df9f5SJed Brown 
1417ab1df9f5SJed Brown     Collective on PC
1418ab1df9f5SJed Brown 
1419ab1df9f5SJed Brown     Input Parameters:
1420ab1df9f5SJed Brown +   pc  - the preconditioner context
1421c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1422ab1df9f5SJed Brown 
1423ab1df9f5SJed Brown     Options Database:
1424c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1425ab1df9f5SJed Brown 
1426ab1df9f5SJed Brown 
1427ab1df9f5SJed Brown     Level: intermediate
1428ab1df9f5SJed Brown 
1429ab1df9f5SJed Brown     Notes:
1430ab1df9f5SJed Brown     The FULL factorization is
1431ab1df9f5SJed Brown 
1432ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1433ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1434ab1df9f5SJed Brown 
14356be4592eSBarry 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,
14366be4592eSBarry 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).
1437ab1df9f5SJed Brown 
14386be4592eSBarry 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
14396be4592eSBarry 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
14406be4592eSBarry 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,
14416be4592eSBarry 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.
1442ab1df9f5SJed Brown 
14436be4592eSBarry 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
14446be4592eSBarry 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).
1445ab1df9f5SJed Brown 
1446ab1df9f5SJed Brown     References:
1447ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1448ab1df9f5SJed Brown 
1449ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1450ab1df9f5SJed Brown 
1451ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1452ab1df9f5SJed Brown @*/
1453c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1454ab1df9f5SJed Brown {
1455ab1df9f5SJed Brown   PetscErrorCode ierr;
1456ab1df9f5SJed Brown 
1457ab1df9f5SJed Brown   PetscFunctionBegin;
1458ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1459c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1460ab1df9f5SJed Brown   PetscFunctionReturn(0);
1461ab1df9f5SJed Brown }
1462ab1df9f5SJed Brown 
1463ab1df9f5SJed Brown #undef __FUNCT__
1464c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1465c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1466ab1df9f5SJed Brown {
1467ab1df9f5SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1468ab1df9f5SJed Brown 
1469ab1df9f5SJed Brown   PetscFunctionBegin;
1470ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1471ab1df9f5SJed Brown   PetscFunctionReturn(0);
1472ab1df9f5SJed Brown }
1473ab1df9f5SJed Brown 
1474ab1df9f5SJed Brown #undef __FUNCT__
147530ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
147630ad9308SMatthew Knepley /*@C
14778c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
147830ad9308SMatthew Knepley 
147930ad9308SMatthew Knepley    Collective on KSP
148030ad9308SMatthew Knepley 
148130ad9308SMatthew Knepley    Input Parameter:
148230ad9308SMatthew Knepley .  pc - the preconditioner context
148330ad9308SMatthew Knepley 
148430ad9308SMatthew Knepley    Output Parameters:
1485a04f6461SBarry Smith +  A00 - the (0,0) block
1486a04f6461SBarry Smith .  A01 - the (0,1) block
1487a04f6461SBarry Smith .  A10 - the (1,0) block
1488a04f6461SBarry Smith -  A11 - the (1,1) block
148930ad9308SMatthew Knepley 
149030ad9308SMatthew Knepley    Level: advanced
149130ad9308SMatthew Knepley 
149230ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
149330ad9308SMatthew Knepley @*/
1494a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
149530ad9308SMatthew Knepley {
149630ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
149730ad9308SMatthew Knepley 
149830ad9308SMatthew Knepley   PetscFunctionBegin;
14990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1500c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1501a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1502a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1503a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1504a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
150530ad9308SMatthew Knepley   PetscFunctionReturn(0);
150630ad9308SMatthew Knepley }
150730ad9308SMatthew Knepley 
150879416396SBarry Smith EXTERN_C_BEGIN
150979416396SBarry Smith #undef __FUNCT__
151079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
15117087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
151279416396SBarry Smith {
151379416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1514e69d4d44SBarry Smith   PetscErrorCode ierr;
151579416396SBarry Smith 
151679416396SBarry Smith   PetscFunctionBegin;
151779416396SBarry Smith   jac->type = type;
15183b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
15193b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
15203b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1521e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1522e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1523c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1524e69d4d44SBarry Smith 
15253b224e63SBarry Smith   } else {
15263b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
15273b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1528e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
15299e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1530c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
15313b224e63SBarry Smith   }
153279416396SBarry Smith   PetscFunctionReturn(0);
153379416396SBarry Smith }
153479416396SBarry Smith EXTERN_C_END
153579416396SBarry Smith 
153651f519a2SBarry Smith EXTERN_C_BEGIN
153751f519a2SBarry Smith #undef __FUNCT__
153851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
15397087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
154051f519a2SBarry Smith {
154151f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
154251f519a2SBarry Smith 
154351f519a2SBarry Smith   PetscFunctionBegin;
154465e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
154565e19b50SBarry 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);
154651f519a2SBarry Smith   jac->bs = bs;
154751f519a2SBarry Smith   PetscFunctionReturn(0);
154851f519a2SBarry Smith }
154951f519a2SBarry Smith EXTERN_C_END
155051f519a2SBarry Smith 
155179416396SBarry Smith #undef __FUNCT__
155279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1553bc08b0f1SBarry Smith /*@
155479416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
155579416396SBarry Smith 
155679416396SBarry Smith    Collective on PC
155779416396SBarry Smith 
155879416396SBarry Smith    Input Parameter:
155979416396SBarry Smith .  pc - the preconditioner context
156081540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
156179416396SBarry Smith 
156279416396SBarry Smith    Options Database Key:
1563a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
156479416396SBarry Smith 
1565b02e2d75SMatthew G Knepley    Level: Intermediate
156679416396SBarry Smith 
156779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
156879416396SBarry Smith 
156979416396SBarry Smith .seealso: PCCompositeSetType()
157079416396SBarry Smith 
157179416396SBarry Smith @*/
15727087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
157379416396SBarry Smith {
15744ac538c5SBarry Smith   PetscErrorCode ierr;
157579416396SBarry Smith 
157679416396SBarry Smith   PetscFunctionBegin;
15770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15784ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
157979416396SBarry Smith   PetscFunctionReturn(0);
158079416396SBarry Smith }
158179416396SBarry Smith 
1582b02e2d75SMatthew G Knepley #undef __FUNCT__
1583b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType"
1584b02e2d75SMatthew G Knepley /*@
1585b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
1586b02e2d75SMatthew G Knepley 
1587b02e2d75SMatthew G Knepley   Not collective
1588b02e2d75SMatthew G Knepley 
1589b02e2d75SMatthew G Knepley   Input Parameter:
1590b02e2d75SMatthew G Knepley . pc - the preconditioner context
1591b02e2d75SMatthew G Knepley 
1592b02e2d75SMatthew G Knepley   Output Parameter:
1593b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1594b02e2d75SMatthew G Knepley 
1595b02e2d75SMatthew G Knepley   Level: Intermediate
1596b02e2d75SMatthew G Knepley 
1597b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1598b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
1599b02e2d75SMatthew G Knepley @*/
1600b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
1601b02e2d75SMatthew G Knepley {
1602b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
1603b02e2d75SMatthew G Knepley 
1604b02e2d75SMatthew G Knepley   PetscFunctionBegin;
1605b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1606b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
1607b02e2d75SMatthew G Knepley   *type = jac->type;
1608b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
1609b02e2d75SMatthew G Knepley }
1610b02e2d75SMatthew G Knepley 
16110971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
16120971522cSBarry Smith /*MC
1613a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1614a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
16150971522cSBarry Smith 
1616edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1617edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
16180971522cSBarry Smith 
1619a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
162069a612a9SBarry Smith          and set the options directly on the resulting KSP object
16210971522cSBarry Smith 
16220971522cSBarry Smith    Level: intermediate
16230971522cSBarry Smith 
162479416396SBarry Smith    Options Database Keys:
162581540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
162681540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
162781540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
162881540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
16290f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
16300f188ba9SJed Brown .   -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag
1631435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
163279416396SBarry Smith 
16335d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
16345d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
16355d4c12cdSJungho Lee 
1636c8a0d604SMatthew G Knepley    Notes:
1637c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1638d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1639d32f9abdSBarry Smith 
1640d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1641d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1642d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1643d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1644d32f9abdSBarry Smith 
1645c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1646c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1647c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1648c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1649c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1650a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1651c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1652c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1653c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1654a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1655c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1656c9c6ffaaSJed Brown      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
16575668aaf4SBarry Smith      diag gives
1658c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1659c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
16605668aaf4SBarry 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
1661c8a0d604SMatthew G Knepley $              (  A00   0 )
1662c8a0d604SMatthew G Knepley $              (  A10   S )
1663c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1664c8a0d604SMatthew G Knepley $              ( A00 A01 )
1665c8a0d604SMatthew G Knepley $              (  0   S  )
1666c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1667e69d4d44SBarry Smith 
1668edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1669edf189efSBarry Smith      is used automatically for a second block.
1670edf189efSBarry Smith 
1671ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1672ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1673ff218e97SBarry Smith 
1674ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1675ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1676ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
16770716a85fSBarry Smith 
1678a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
16790971522cSBarry Smith 
16807e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1681e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
16820971522cSBarry Smith M*/
16830971522cSBarry Smith 
16840971522cSBarry Smith EXTERN_C_BEGIN
16850971522cSBarry Smith #undef __FUNCT__
16860971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
16877087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
16880971522cSBarry Smith {
16890971522cSBarry Smith   PetscErrorCode ierr;
16900971522cSBarry Smith   PC_FieldSplit  *jac;
16910971522cSBarry Smith 
16920971522cSBarry Smith   PetscFunctionBegin;
169338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
16940971522cSBarry Smith   jac->bs        = -1;
16950971522cSBarry Smith   jac->nsplits   = 0;
16963e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1697e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1698c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
169951f519a2SBarry Smith 
17000971522cSBarry Smith   pc->data     = (void*)jac;
17010971522cSBarry Smith 
17020971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1703421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
17040971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1705574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
17060971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
17070971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
17080971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
17090971522cSBarry Smith   pc->ops->applyrichardson   = 0;
17100971522cSBarry Smith 
171169a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
171269a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
17130971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
17140971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1715704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1716704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
171779416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
171879416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
171951f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
172051f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
17210971522cSBarry Smith   PetscFunctionReturn(0);
17220971522cSBarry Smith }
17230971522cSBarry Smith EXTERN_C_END
17240971522cSBarry Smith 
17250971522cSBarry Smith 
1726a541d17aSBarry Smith 
1727