xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1dba47a55SKris Buschelman 
2*b45d2f2cSJed 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;
762692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
770971522cSBarry Smith   if (iascii) {
7851f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
7969a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
800971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
810971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
821ab39975SBarry Smith       if (ilink->fields) {
830971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
8479416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
855a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
8679416396SBarry Smith 	  if (j > 0) {
8779416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
8879416396SBarry Smith 	  }
895a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
900971522cSBarry Smith 	}
910971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9279416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
931ab39975SBarry Smith       } else {
941ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
951ab39975SBarry Smith       }
965a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
975a9f2f41SSatish Balay       ilink = ilink->next;
980971522cSBarry Smith     }
990971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1000971522cSBarry Smith   } else {
10165e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1020971522cSBarry Smith   }
1030971522cSBarry Smith   PetscFunctionReturn(0);
1040971522cSBarry Smith }
1050971522cSBarry Smith 
1060971522cSBarry Smith #undef __FUNCT__
1073b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1083b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1093b224e63SBarry Smith {
1103b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1113b224e63SBarry Smith   PetscErrorCode    ierr;
112ace3abfcSBarry Smith   PetscBool         iascii;
1133b224e63SBarry Smith   PetscInt          i,j;
1143b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1153b224e63SBarry Smith   KSP               ksp;
1163b224e63SBarry Smith 
1173b224e63SBarry Smith   PetscFunctionBegin;
1182692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1193b224e63SBarry Smith   if (iascii) {
120c9c6ffaaSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1213b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1223b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1233b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1243b224e63SBarry Smith       if (ilink->fields) {
1253b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1263b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1273b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1283b224e63SBarry Smith 	  if (j > 0) {
1293b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1303b224e63SBarry Smith 	  }
1313b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1323b224e63SBarry Smith 	}
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1343b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1353b224e63SBarry Smith       } else {
1363b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1373b224e63SBarry Smith       }
1383b224e63SBarry Smith       ilink = ilink->next;
1393b224e63SBarry Smith     }
140435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1413b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14212cae6f2SJed Brown     if (jac->schur) {
14312cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1443b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
14512cae6f2SJed Brown     } else {
14612cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
14712cae6f2SJed Brown     }
1483b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
149435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1503b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15112cae6f2SJed Brown     if (jac->kspschur) {
1523b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
15312cae6f2SJed Brown     } else {
15412cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15512cae6f2SJed Brown     }
1563b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1573b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1583b224e63SBarry Smith   } else {
15965e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1603b224e63SBarry Smith   }
1613b224e63SBarry Smith   PetscFunctionReturn(0);
1623b224e63SBarry Smith }
1633b224e63SBarry Smith 
1643b224e63SBarry Smith #undef __FUNCT__
1656c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1666c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1676c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1686c924f48SJed Brown {
1696c924f48SJed Brown   PetscErrorCode ierr;
1706c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1715d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
1725d4c12cdSJungho Lee   PetscBool      flg,flg_col;
1735d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
1746c924f48SJed Brown 
1756c924f48SJed Brown   PetscFunctionBegin;
1766c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1775d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
1786c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1796c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1806c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1815d4c12cdSJungho Lee     ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
1826c924f48SJed Brown     nfields = jac->bs;
18329499fbbSJungho Lee     nfields_col = jac->bs;
1846c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1855d4c12cdSJungho Lee     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
1866c924f48SJed Brown     if (!flg) break;
1875d4c12cdSJungho Lee     else if (flg && !flg_col) {
1885d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1895d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
1905d4c12cdSJungho Lee     }
1915d4c12cdSJungho Lee     else {
1925d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1935d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
1945d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
1955d4c12cdSJungho Lee     }
1966c924f48SJed Brown   }
1976c924f48SJed Brown   if (i > 0) {
1986c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1996c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2006c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2016c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2026c924f48SJed Brown   }
2036c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2045d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2056c924f48SJed Brown   PetscFunctionReturn(0);
2066c924f48SJed Brown }
2076c924f48SJed Brown 
2086c924f48SJed Brown #undef __FUNCT__
20969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
21069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2110971522cSBarry Smith {
2120971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2130971522cSBarry Smith   PetscErrorCode    ierr;
2145a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2156ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2166c924f48SJed Brown   PetscInt          i;
2170971522cSBarry Smith 
2180971522cSBarry Smith   PetscFunctionBegin;
219d32f9abdSBarry Smith   if (!ilink) {
220d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
221d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2227b62db95SJungho Lee       PetscInt     numFields, f;
2230784a22cSJed Brown       char         **fieldNames;
2247b62db95SJungho Lee       IS          *fields;
2254d343eeaSMatthew G Knepley       PetscBool    dmcomposite;
2267b62db95SJungho Lee 
2277b62db95SJungho Lee       ierr = DMCreateFieldIS(pc->dm, &numFields, &fieldNames, &fields);CHKERRQ(ierr);
2287b62db95SJungho Lee       for(f = 0; f < numFields; ++f) {
2297b62db95SJungho Lee         ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
2307b62db95SJungho Lee         ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
2317b62db95SJungho Lee         ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
2327b62db95SJungho Lee       }
2337b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
2347b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
2357b62db95SJungho Lee 
2364d343eeaSMatthew G Knepley       ierr = PetscTypeCompare((PetscObject) pc->dm, DMCOMPOSITE, &dmcomposite);CHKERRQ(ierr);
2374d343eeaSMatthew G Knepley       if (dmcomposite) {
2387b62db95SJungho Lee         DM      *dms;
2394d343eeaSMatthew G Knepley         PetscInt nDM;
2407b62db95SJungho Lee 
2418b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2428b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm, &nDM);CHKERRQ(ierr);
2437b62db95SJungho Lee         ierr = PetscMalloc(nDM*sizeof(DM), &dms);CHKERRQ(ierr);
2447b62db95SJungho Lee         ierr = DMCompositeGetEntriesArray(pc->dm, dms);CHKERRQ(ierr);
2457b62db95SJungho Lee         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2467b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
2477b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
2482fa5ba8aSJed Brown         }
2497b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
2508b8307b2SJed Brown       }
25166ffff09SJed Brown     } else {
252521d7252SBarry Smith       if (jac->bs <= 0) {
253704ba839SBarry Smith         if (pc->pmat) {
254521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
255704ba839SBarry Smith         } else {
256704ba839SBarry Smith           jac->bs = 1;
257704ba839SBarry Smith         }
258521d7252SBarry Smith       }
259d32f9abdSBarry Smith 
260acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2616ce1633cSBarry Smith       if (stokes) {
2626ce1633cSBarry Smith         IS       zerodiags,rest;
2636ce1633cSBarry Smith         PetscInt nmin,nmax;
2646ce1633cSBarry Smith 
2656ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2666ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2676ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2686ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2696ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
270fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
271fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2726ce1633cSBarry Smith       } else {
273d32f9abdSBarry Smith         if (!flg) {
274d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
275d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2766c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2776c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
278d32f9abdSBarry Smith         }
2796c924f48SJed Brown         if (flg || !jac->splitdefined) {
280d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
281db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2826c924f48SJed Brown             char splitname[8];
2836c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2845d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
28579416396SBarry Smith           }
2865d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
287521d7252SBarry Smith         }
28866ffff09SJed Brown       }
2896ce1633cSBarry Smith     }
290edf189efSBarry Smith   } else if (jac->nsplits == 1) {
291edf189efSBarry Smith     if (ilink->is) {
292edf189efSBarry Smith       IS       is2;
293edf189efSBarry Smith       PetscInt nmin,nmax;
294edf189efSBarry Smith 
295edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
296edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
297db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
298fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
2997b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
30063ec74ffSBarry Smith   } else if (jac->reset) {
30163ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
302d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
303d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
304d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
305d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
306d0af7cd3SBarry Smith       PetscBool dmcomposite;
307d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
308d0af7cd3SBarry Smith       if (dmcomposite) {
309d0af7cd3SBarry Smith         PetscInt nDM;
310d0af7cd3SBarry Smith         IS       *fields;
3117b62db95SJungho Lee         DM       *dms;
312d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
313d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
314d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3157b62db95SJungho Lee         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3167b62db95SJungho Lee         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
317d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3187b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3197b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
320d0af7cd3SBarry Smith           ilink->is = fields[i];
321d0af7cd3SBarry Smith           ilink     = ilink->next;
322edf189efSBarry Smith         }
323d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3247b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
325d0af7cd3SBarry Smith       }
326d0af7cd3SBarry Smith     } else {
327d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
328d0af7cd3SBarry Smith       if (stokes) {
329d0af7cd3SBarry Smith         IS       zerodiags,rest;
330d0af7cd3SBarry Smith         PetscInt nmin,nmax;
331d0af7cd3SBarry Smith 
332d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
333d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
334d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
33520b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
33620b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
337d0af7cd3SBarry Smith         ilink->is       = rest;
338d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
33920b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
340d0af7cd3SBarry Smith     }
341d0af7cd3SBarry Smith   }
342d0af7cd3SBarry Smith 
3437b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
34469a612a9SBarry Smith   PetscFunctionReturn(0);
34569a612a9SBarry Smith }
34669a612a9SBarry Smith 
34769a612a9SBarry Smith #undef __FUNCT__
34869a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
34969a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
35069a612a9SBarry Smith {
35169a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
35269a612a9SBarry Smith   PetscErrorCode    ierr;
3535a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
354cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
35569a612a9SBarry Smith   MatStructure      flag = pc->flag;
3565f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
35769a612a9SBarry Smith 
35869a612a9SBarry Smith   PetscFunctionBegin;
35969a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
36097bbdb24SBarry Smith   nsplit = jac->nsplits;
3615a9f2f41SSatish Balay   ilink  = jac->head;
36297bbdb24SBarry Smith 
36397bbdb24SBarry Smith   /* get the matrices for each split */
364704ba839SBarry Smith   if (!jac->issetup) {
3651b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
36697bbdb24SBarry Smith 
367704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
368704ba839SBarry Smith 
3695d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
37051f519a2SBarry Smith     bs     = jac->bs;
37197bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
372cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3731b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3741b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3751b9fc7fcSBarry Smith       if (jac->defaultsplit) {
376704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
3775f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
378704ba839SBarry Smith       } else if (!ilink->is) {
379ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3805f4ab4e1SJungho Lee           PetscInt   *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
3815a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3825f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
3831b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3841b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3851b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
3865f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
38797bbdb24SBarry Smith             }
38897bbdb24SBarry Smith           }
389d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
3905f4ab4e1SJungho Lee           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
391ccb205f8SBarry Smith         } else {
392704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
3935f4ab4e1SJungho Lee           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
394ccb205f8SBarry Smith        }
3953e197d65SBarry Smith       }
396704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
3975f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
3985f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
399704ba839SBarry Smith       ilink = ilink->next;
4001b9fc7fcSBarry Smith     }
4011b9fc7fcSBarry Smith   }
4021b9fc7fcSBarry Smith 
403704ba839SBarry Smith   ilink  = jac->head;
40497bbdb24SBarry Smith   if (!jac->pmat) {
405cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
406cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4075f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
408704ba839SBarry Smith       ilink = ilink->next;
409cf502942SBarry Smith     }
41097bbdb24SBarry Smith   } else {
411cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4125f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
413704ba839SBarry Smith       ilink = ilink->next;
414cf502942SBarry Smith     }
41597bbdb24SBarry Smith   }
416519d70e2SJed Brown   if (jac->realdiagonal) {
417519d70e2SJed Brown     ilink = jac->head;
418519d70e2SJed Brown     if (!jac->mat) {
419519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
420519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4215f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
422519d70e2SJed Brown         ilink = ilink->next;
423519d70e2SJed Brown       }
424519d70e2SJed Brown     } else {
425519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4265f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
427519d70e2SJed Brown         ilink = ilink->next;
428519d70e2SJed Brown       }
429519d70e2SJed Brown     }
430519d70e2SJed Brown   } else {
431519d70e2SJed Brown     jac->mat = jac->pmat;
432519d70e2SJed Brown   }
43397bbdb24SBarry Smith 
4346c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
43568dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
43668dd23aaSBarry Smith     ilink  = jac->head;
43768dd23aaSBarry Smith     if (!jac->Afield) {
43868dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
43968dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4404aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
44168dd23aaSBarry Smith         ilink = ilink->next;
44268dd23aaSBarry Smith       }
44368dd23aaSBarry Smith     } else {
44468dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4454aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
44668dd23aaSBarry Smith         ilink = ilink->next;
44768dd23aaSBarry Smith       }
44868dd23aaSBarry Smith     }
44968dd23aaSBarry Smith   }
45068dd23aaSBarry Smith 
4513b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4523b224e63SBarry Smith     IS       ccis;
4534aa3045dSJed Brown     PetscInt rstart,rend;
454093c86ffSJed Brown     char     lscname[256];
455093c86ffSJed Brown     PetscObject LSC_L;
456e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
45768dd23aaSBarry Smith 
458e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
459e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
460e6cab6aaSJed Brown 
4613b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4623b224e63SBarry Smith     if (jac->schur) {
4633b224e63SBarry Smith       ilink = jac->head;
46449bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4654aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
466fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4673b224e63SBarry Smith       ilink = ilink->next;
46849bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4694aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
470fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
471519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
472084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4733b224e63SBarry Smith 
4743b224e63SBarry Smith      } else {
4751cee3971SBarry Smith       KSP ksp;
4766c924f48SJed Brown       char schurprefix[256];
4773b224e63SBarry Smith 
478a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4793b224e63SBarry Smith       ilink = jac->head;
48049bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4814aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
482fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4833b224e63SBarry Smith       ilink = ilink->next;
48449bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4854aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
486fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4877d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
488519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
489a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4901cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
491e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
492a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
493a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
49420b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
49520b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4967b62db95SJungho Lee       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
4973b224e63SBarry Smith 
4983b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4999005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
5001cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
501084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
502084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
503084e4875SJed Brown         PC pc;
504cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
505084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
506084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
507e69d4d44SBarry Smith       }
5086c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5096c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5103b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
51120b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
51220b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5133b224e63SBarry Smith 
5143b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5153b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5163b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5173b224e63SBarry Smith       ilink = jac->head;
5183b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5193b224e63SBarry Smith       ilink = ilink->next;
5203b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5213b224e63SBarry Smith     }
522093c86ffSJed Brown 
523093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
524093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
525093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
526093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
527093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
528093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
529093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
530093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
531093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
5323b224e63SBarry Smith   } else {
53397bbdb24SBarry Smith     /* set up the individual PCs */
53497bbdb24SBarry Smith     i    = 0;
5355a9f2f41SSatish Balay     ilink = jac->head;
5365a9f2f41SSatish Balay     while (ilink) {
537519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5383b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
539c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5405a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
54197bbdb24SBarry Smith       i++;
5425a9f2f41SSatish Balay       ilink = ilink->next;
5430971522cSBarry Smith     }
54497bbdb24SBarry Smith 
54597bbdb24SBarry Smith     /* create work vectors for each split */
5461b9fc7fcSBarry Smith     if (!jac->x) {
54797bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5485a9f2f41SSatish Balay       ilink = jac->head;
54997bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
550906ed7ccSBarry Smith         Vec *vl,*vr;
551906ed7ccSBarry Smith 
552906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
553906ed7ccSBarry Smith         ilink->x  = *vr;
554906ed7ccSBarry Smith         ilink->y  = *vl;
555906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
556906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5575a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5585a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5595a9f2f41SSatish Balay         ilink     = ilink->next;
56097bbdb24SBarry Smith       }
5613b224e63SBarry Smith     }
5623b224e63SBarry Smith   }
5633b224e63SBarry Smith 
5643b224e63SBarry Smith 
565d0f46423SBarry Smith   if (!jac->head->sctx) {
5663b224e63SBarry Smith     Vec xtmp;
5673b224e63SBarry Smith 
56879416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5691b9fc7fcSBarry Smith 
5705a9f2f41SSatish Balay     ilink = jac->head;
5711b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5721b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
573704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5745a9f2f41SSatish Balay       ilink = ilink->next;
5751b9fc7fcSBarry Smith     }
5766bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5771b9fc7fcSBarry Smith   }
578c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5790971522cSBarry Smith   PetscFunctionReturn(0);
5800971522cSBarry Smith }
5810971522cSBarry Smith 
5825a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
583ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
584ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5855a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
586ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
587ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
58879416396SBarry Smith 
5890971522cSBarry Smith #undef __FUNCT__
5903b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5913b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5923b224e63SBarry Smith {
5933b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5943b224e63SBarry Smith   PetscErrorCode    ierr;
5953b224e63SBarry Smith   KSP               ksp;
5963b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5973b224e63SBarry Smith 
5983b224e63SBarry Smith   PetscFunctionBegin;
5993b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
6003b224e63SBarry Smith 
601c5d2311dSJed Brown   switch (jac->schurfactorization) {
602c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
603a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
604c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
605c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
606c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
607c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
608c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
609c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
610c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
611c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
615c5d2311dSJed Brown       break;
616c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
617a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
618c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
619c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
620c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
621c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
622c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
623c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
624c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
625c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
627c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
628c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
629c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
630c5d2311dSJed Brown       break;
631c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
632a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
633c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
634c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
635c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
636c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
637c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
638c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
639c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
640c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
641c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
642c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
643c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
644c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
645c5d2311dSJed Brown       break;
646c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_FULL:
647a04f6461SBarry 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 */
6483b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6493b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6503b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6513b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6523b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6533b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6543b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6553b224e63SBarry Smith 
6563b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6573b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6583b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6593b224e63SBarry Smith 
6603b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6613b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6623b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6633b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6643b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
665c5d2311dSJed Brown   }
6663b224e63SBarry Smith   PetscFunctionReturn(0);
6673b224e63SBarry Smith }
6683b224e63SBarry Smith 
6693b224e63SBarry Smith #undef __FUNCT__
6700971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6710971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6720971522cSBarry Smith {
6730971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6740971522cSBarry Smith   PetscErrorCode    ierr;
6755a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
676af93645bSJed Brown   PetscInt          cnt;
6770971522cSBarry Smith 
6780971522cSBarry Smith   PetscFunctionBegin;
67951f519a2SBarry Smith   CHKMEMQ;
68051f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
68151f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
68251f519a2SBarry Smith 
68379416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6841b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6850971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6865a9f2f41SSatish Balay       while (ilink) {
6875a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6885a9f2f41SSatish Balay         ilink = ilink->next;
6890971522cSBarry Smith       }
6900971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6911b9fc7fcSBarry Smith     } else {
692efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6935a9f2f41SSatish Balay       while (ilink) {
6945a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6955a9f2f41SSatish Balay         ilink = ilink->next;
6961b9fc7fcSBarry Smith       }
6971b9fc7fcSBarry Smith     }
69816913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
69979416396SBarry Smith     if (!jac->w1) {
70079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
70179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
70279416396SBarry Smith     }
703efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
7045a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7053e197d65SBarry Smith     cnt = 1;
7065a9f2f41SSatish Balay     while (ilink->next) {
7075a9f2f41SSatish Balay       ilink = ilink->next;
7083e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7093e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7103e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7113e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7123e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7133e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7143e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7153e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7163e197d65SBarry Smith     }
71751f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
71811755939SBarry Smith       cnt -= 2;
71951f519a2SBarry Smith       while (ilink->previous) {
72051f519a2SBarry Smith         ilink = ilink->previous;
72111755939SBarry Smith         /* compute the residual only over the part of the vector needed */
72211755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
72311755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
72411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72611755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
72711755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72811755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72951f519a2SBarry Smith       }
73011755939SBarry Smith     }
73165e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
73251f519a2SBarry Smith   CHKMEMQ;
7330971522cSBarry Smith   PetscFunctionReturn(0);
7340971522cSBarry Smith }
7350971522cSBarry Smith 
736421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
737ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
738ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
739421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
740ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
741ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
742421e10b8SBarry Smith 
743421e10b8SBarry Smith #undef __FUNCT__
7448c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
745421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
746421e10b8SBarry Smith {
747421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
748421e10b8SBarry Smith   PetscErrorCode    ierr;
749421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
750421e10b8SBarry Smith 
751421e10b8SBarry Smith   PetscFunctionBegin;
752421e10b8SBarry Smith   CHKMEMQ;
753421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
754421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
755421e10b8SBarry Smith 
756421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
757421e10b8SBarry Smith     if (jac->defaultsplit) {
758421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
759421e10b8SBarry Smith       while (ilink) {
760421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
761421e10b8SBarry Smith 	ilink = ilink->next;
762421e10b8SBarry Smith       }
763421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
764421e10b8SBarry Smith     } else {
765421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
766421e10b8SBarry Smith       while (ilink) {
767421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
768421e10b8SBarry Smith 	ilink = ilink->next;
769421e10b8SBarry Smith       }
770421e10b8SBarry Smith     }
771421e10b8SBarry Smith   } else {
772421e10b8SBarry Smith     if (!jac->w1) {
773421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
774421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
775421e10b8SBarry Smith     }
776421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
777421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
778421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
779421e10b8SBarry Smith       while (ilink->next) {
780421e10b8SBarry Smith         ilink = ilink->next;
7819989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
782421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
783421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
784421e10b8SBarry Smith       }
785421e10b8SBarry Smith       while (ilink->previous) {
786421e10b8SBarry Smith         ilink = ilink->previous;
7879989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
788421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
789421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
790421e10b8SBarry Smith       }
791421e10b8SBarry Smith     } else {
792421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
793421e10b8SBarry Smith 	ilink = ilink->next;
794421e10b8SBarry Smith       }
795421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
796421e10b8SBarry Smith       while (ilink->previous) {
797421e10b8SBarry Smith 	ilink = ilink->previous;
7989989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
799421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
800421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
801421e10b8SBarry Smith       }
802421e10b8SBarry Smith     }
803421e10b8SBarry Smith   }
804421e10b8SBarry Smith   CHKMEMQ;
805421e10b8SBarry Smith   PetscFunctionReturn(0);
806421e10b8SBarry Smith }
807421e10b8SBarry Smith 
8080971522cSBarry Smith #undef __FUNCT__
809574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
810574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8110971522cSBarry Smith {
8120971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8130971522cSBarry Smith   PetscErrorCode    ierr;
8145a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8150971522cSBarry Smith 
8160971522cSBarry Smith   PetscFunctionBegin;
8175a9f2f41SSatish Balay   while (ilink) {
818574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
819fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
820fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
821fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
822fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
823b5787286SJed Brown     ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
8245a9f2f41SSatish Balay     next = ilink->next;
8255a9f2f41SSatish Balay     ilink = next;
8260971522cSBarry Smith   }
82705b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
828574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
829574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
830574deadeSBarry Smith   } else if (jac->mat) {
831574deadeSBarry Smith     jac->mat = PETSC_NULL;
832574deadeSBarry Smith   }
83397bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
83468dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8356bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8366bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8376bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8386bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8396bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8406bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8416bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
84263ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
843574deadeSBarry Smith   PetscFunctionReturn(0);
844574deadeSBarry Smith }
845574deadeSBarry Smith 
846574deadeSBarry Smith #undef __FUNCT__
847574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
848574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
849574deadeSBarry Smith {
850574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
851574deadeSBarry Smith   PetscErrorCode    ierr;
852574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
853574deadeSBarry Smith 
854574deadeSBarry Smith   PetscFunctionBegin;
855574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
856574deadeSBarry Smith   while (ilink) {
8576bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
858574deadeSBarry Smith     next = ilink->next;
859574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
860574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
8615d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
862574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
863574deadeSBarry Smith     ilink = next;
864574deadeSBarry Smith   }
865574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
866c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8677b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
8687b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
8697b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
8707b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
8717b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
872ab1df9f5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr);
873c9c6ffaaSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr);
8740971522cSBarry Smith   PetscFunctionReturn(0);
8750971522cSBarry Smith }
8760971522cSBarry Smith 
8770971522cSBarry Smith #undef __FUNCT__
8780971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8790971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8800971522cSBarry Smith {
8811b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8826c924f48SJed Brown   PetscInt        bs;
883bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8849dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8853b224e63SBarry Smith   PCCompositeType ctype;
8861b9fc7fcSBarry Smith 
8870971522cSBarry Smith   PetscFunctionBegin;
8881b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
889acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
89051f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
89151f519a2SBarry Smith   if (flg) {
89251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
89351f519a2SBarry Smith   }
894704ba839SBarry Smith 
895435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
896c0adfefeSBarry Smith   if (stokes) {
897c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
898c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
899c0adfefeSBarry Smith   }
900c0adfefeSBarry Smith 
9013b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
9023b224e63SBarry Smith   if (flg) {
9033b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
9043b224e63SBarry Smith   }
905d32f9abdSBarry Smith 
906c30613efSMatthew Knepley   /* Only setup fields once */
907c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
908d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
909d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
9106c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9116c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
912d32f9abdSBarry Smith   }
913c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
914c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
915c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
916c9c6ffaaSJed 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);
917084e4875SJed 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);
918c5d2311dSJed Brown   }
9191b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9200971522cSBarry Smith   PetscFunctionReturn(0);
9210971522cSBarry Smith }
9220971522cSBarry Smith 
9230971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9240971522cSBarry Smith 
9250971522cSBarry Smith EXTERN_C_BEGIN
9260971522cSBarry Smith #undef __FUNCT__
9270971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9285d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
9290971522cSBarry Smith {
93097bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9310971522cSBarry Smith   PetscErrorCode    ierr;
9325a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
93369a612a9SBarry Smith   char              prefix[128];
9345d4c12cdSJungho Lee   PetscInt          i;
9350971522cSBarry Smith 
9360971522cSBarry Smith   PetscFunctionBegin;
9376c924f48SJed Brown   if (jac->splitdefined) {
9386c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9396c924f48SJed Brown     PetscFunctionReturn(0);
9406c924f48SJed Brown   }
94151f519a2SBarry Smith   for (i=0; i<n; i++) {
942e32f2f54SBarry 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);
943e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
94451f519a2SBarry Smith   }
945704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
946a04f6461SBarry Smith   if (splitname) {
947db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
948a04f6461SBarry Smith   } else {
949a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
950a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
951a04f6461SBarry Smith   }
952704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9535a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9545d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
9555d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
9565a9f2f41SSatish Balay   ilink->nfields = n;
9575a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9587adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9591cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9605a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9619005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
96269a612a9SBarry Smith 
963a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9645a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9650971522cSBarry Smith 
9660971522cSBarry Smith   if (!next) {
9675a9f2f41SSatish Balay     jac->head       = ilink;
96851f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9690971522cSBarry Smith   } else {
9700971522cSBarry Smith     while (next->next) {
9710971522cSBarry Smith       next = next->next;
9720971522cSBarry Smith     }
9735a9f2f41SSatish Balay     next->next      = ilink;
97451f519a2SBarry Smith     ilink->previous = next;
9750971522cSBarry Smith   }
9760971522cSBarry Smith   jac->nsplits++;
9770971522cSBarry Smith   PetscFunctionReturn(0);
9780971522cSBarry Smith }
9790971522cSBarry Smith EXTERN_C_END
9800971522cSBarry Smith 
981e69d4d44SBarry Smith EXTERN_C_BEGIN
982e69d4d44SBarry Smith #undef __FUNCT__
983e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9847087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
985e69d4d44SBarry Smith {
986e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
987e69d4d44SBarry Smith   PetscErrorCode ierr;
988e69d4d44SBarry Smith 
989e69d4d44SBarry Smith   PetscFunctionBegin;
990519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
991e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
992e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
99313e0d083SBarry Smith   if (n) *n = jac->nsplits;
994e69d4d44SBarry Smith   PetscFunctionReturn(0);
995e69d4d44SBarry Smith }
996e69d4d44SBarry Smith EXTERN_C_END
9970971522cSBarry Smith 
9980971522cSBarry Smith EXTERN_C_BEGIN
9990971522cSBarry Smith #undef __FUNCT__
100069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
10017087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
10020971522cSBarry Smith {
10030971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10040971522cSBarry Smith   PetscErrorCode    ierr;
10050971522cSBarry Smith   PetscInt          cnt = 0;
10065a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
10070971522cSBarry Smith 
10080971522cSBarry Smith   PetscFunctionBegin;
10095d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
10105a9f2f41SSatish Balay   while (ilink) {
10115a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
10125a9f2f41SSatish Balay     ilink = ilink->next;
10130971522cSBarry Smith   }
10145d480477SMatthew 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);
101513e0d083SBarry Smith   if (n) *n = jac->nsplits;
10160971522cSBarry Smith   PetscFunctionReturn(0);
10170971522cSBarry Smith }
10180971522cSBarry Smith EXTERN_C_END
10190971522cSBarry Smith 
1020704ba839SBarry Smith EXTERN_C_BEGIN
1021704ba839SBarry Smith #undef __FUNCT__
1022704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10237087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1024704ba839SBarry Smith {
1025704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1026704ba839SBarry Smith   PetscErrorCode    ierr;
1027704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1028704ba839SBarry Smith   char              prefix[128];
1029704ba839SBarry Smith 
1030704ba839SBarry Smith   PetscFunctionBegin;
10316c924f48SJed Brown   if (jac->splitdefined) {
10326c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10336c924f48SJed Brown     PetscFunctionReturn(0);
10346c924f48SJed Brown   }
103516913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1036a04f6461SBarry Smith   if (splitname) {
1037db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1038a04f6461SBarry Smith   } else {
1039b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1040b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1041a04f6461SBarry Smith   }
10421ab39975SBarry Smith   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1043b5787286SJed Brown   ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
1044b5787286SJed Brown   ilink->is      = is;
1045b5787286SJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1046b5787286SJed Brown   ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1047b5787286SJed Brown   ilink->is_col  = is;
1048704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1049704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10501cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1051704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10529005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1053704ba839SBarry Smith 
1054a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1055704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1056704ba839SBarry Smith 
1057704ba839SBarry Smith   if (!next) {
1058704ba839SBarry Smith     jac->head       = ilink;
1059704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1060704ba839SBarry Smith   } else {
1061704ba839SBarry Smith     while (next->next) {
1062704ba839SBarry Smith       next = next->next;
1063704ba839SBarry Smith     }
1064704ba839SBarry Smith     next->next      = ilink;
1065704ba839SBarry Smith     ilink->previous = next;
1066704ba839SBarry Smith   }
1067704ba839SBarry Smith   jac->nsplits++;
1068704ba839SBarry Smith 
1069704ba839SBarry Smith   PetscFunctionReturn(0);
1070704ba839SBarry Smith }
1071704ba839SBarry Smith EXTERN_C_END
1072704ba839SBarry Smith 
10730971522cSBarry Smith #undef __FUNCT__
10740971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10750971522cSBarry Smith /*@
10760971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10770971522cSBarry Smith 
1078ad4df100SBarry Smith     Logically Collective on PC
10790971522cSBarry Smith 
10800971522cSBarry Smith     Input Parameters:
10810971522cSBarry Smith +   pc  - the preconditioner context
1082a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10830971522cSBarry Smith .   n - the number of fields in this split
1084db4c96c1SJed Brown -   fields - the fields in this split
10850971522cSBarry Smith 
10860971522cSBarry Smith     Level: intermediate
10870971522cSBarry Smith 
1088d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1089d32f9abdSBarry Smith 
1090d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1091d32f9abdSBarry 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
1092d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1093d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1094d32f9abdSBarry Smith 
1095db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1096db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1097db4c96c1SJed Brown 
10985d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
10995d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
11005d4c12cdSJungho Lee      available when this routine is called.
11015d4c12cdSJungho Lee 
1102d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
11030971522cSBarry Smith 
11040971522cSBarry Smith @*/
11055d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
11060971522cSBarry Smith {
11074ac538c5SBarry Smith   PetscErrorCode ierr;
11080971522cSBarry Smith 
11090971522cSBarry Smith   PetscFunctionBegin;
11100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1111db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1112db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1113db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
11145d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
11150971522cSBarry Smith   PetscFunctionReturn(0);
11160971522cSBarry Smith }
11170971522cSBarry Smith 
11180971522cSBarry Smith #undef __FUNCT__
1119704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1120704ba839SBarry Smith /*@
1121704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1122704ba839SBarry Smith 
1123ad4df100SBarry Smith     Logically Collective on PC
1124704ba839SBarry Smith 
1125704ba839SBarry Smith     Input Parameters:
1126704ba839SBarry Smith +   pc  - the preconditioner context
1127a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1128db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1129704ba839SBarry Smith 
1130d32f9abdSBarry Smith 
1131a6ffb8dbSJed Brown     Notes:
1132a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1133a6ffb8dbSJed Brown 
1134db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1135db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1136d32f9abdSBarry Smith 
1137704ba839SBarry Smith     Level: intermediate
1138704ba839SBarry Smith 
1139704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1140704ba839SBarry Smith 
1141704ba839SBarry Smith @*/
11427087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1143704ba839SBarry Smith {
11444ac538c5SBarry Smith   PetscErrorCode ierr;
1145704ba839SBarry Smith 
1146704ba839SBarry Smith   PetscFunctionBegin;
11470700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11487b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1149db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11504ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1151704ba839SBarry Smith   PetscFunctionReturn(0);
1152704ba839SBarry Smith }
1153704ba839SBarry Smith 
1154704ba839SBarry Smith #undef __FUNCT__
115557a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
115657a9adfeSMatthew G Knepley /*@
115757a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
115857a9adfeSMatthew G Knepley 
115957a9adfeSMatthew G Knepley     Logically Collective on PC
116057a9adfeSMatthew G Knepley 
116157a9adfeSMatthew G Knepley     Input Parameters:
116257a9adfeSMatthew G Knepley +   pc  - the preconditioner context
116357a9adfeSMatthew G Knepley -   splitname - name of this split
116457a9adfeSMatthew G Knepley 
116557a9adfeSMatthew G Knepley     Output Parameter:
116657a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
116757a9adfeSMatthew G Knepley 
116857a9adfeSMatthew G Knepley     Level: intermediate
116957a9adfeSMatthew G Knepley 
117057a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
117157a9adfeSMatthew G Knepley 
117257a9adfeSMatthew G Knepley @*/
117357a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
117457a9adfeSMatthew G Knepley {
117557a9adfeSMatthew G Knepley   PetscErrorCode ierr;
117657a9adfeSMatthew G Knepley 
117757a9adfeSMatthew G Knepley   PetscFunctionBegin;
117857a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
117957a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
118057a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
118157a9adfeSMatthew G Knepley   {
118257a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
118357a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
118457a9adfeSMatthew G Knepley     PetscBool         found;
118557a9adfeSMatthew G Knepley 
118657a9adfeSMatthew G Knepley     *is = PETSC_NULL;
118757a9adfeSMatthew G Knepley     while(ilink) {
118857a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
118957a9adfeSMatthew G Knepley       if (found) {
119057a9adfeSMatthew G Knepley         *is = ilink->is;
119157a9adfeSMatthew G Knepley         break;
119257a9adfeSMatthew G Knepley       }
119357a9adfeSMatthew G Knepley       ilink = ilink->next;
119457a9adfeSMatthew G Knepley     }
119557a9adfeSMatthew G Knepley   }
119657a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
119757a9adfeSMatthew G Knepley }
119857a9adfeSMatthew G Knepley 
119957a9adfeSMatthew G Knepley #undef __FUNCT__
120051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
120151f519a2SBarry Smith /*@
120251f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
120351f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
120451f519a2SBarry Smith 
1205ad4df100SBarry Smith     Logically Collective on PC
120651f519a2SBarry Smith 
120751f519a2SBarry Smith     Input Parameters:
120851f519a2SBarry Smith +   pc  - the preconditioner context
120951f519a2SBarry Smith -   bs - the block size
121051f519a2SBarry Smith 
121151f519a2SBarry Smith     Level: intermediate
121251f519a2SBarry Smith 
121351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
121451f519a2SBarry Smith 
121551f519a2SBarry Smith @*/
12167087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
121751f519a2SBarry Smith {
12184ac538c5SBarry Smith   PetscErrorCode ierr;
121951f519a2SBarry Smith 
122051f519a2SBarry Smith   PetscFunctionBegin;
12210700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1222c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
12234ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
122451f519a2SBarry Smith   PetscFunctionReturn(0);
122551f519a2SBarry Smith }
122651f519a2SBarry Smith 
122751f519a2SBarry Smith #undef __FUNCT__
122869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12290971522cSBarry Smith /*@C
123069a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12310971522cSBarry Smith 
123269a612a9SBarry Smith    Collective on KSP
12330971522cSBarry Smith 
12340971522cSBarry Smith    Input Parameter:
12350971522cSBarry Smith .  pc - the preconditioner context
12360971522cSBarry Smith 
12370971522cSBarry Smith    Output Parameters:
123813e0d083SBarry Smith +  n - the number of splits
123969a612a9SBarry Smith -  pc - the array of KSP contexts
12400971522cSBarry Smith 
12410971522cSBarry Smith    Note:
1242d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1243d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12440971522cSBarry Smith 
124569a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12460971522cSBarry Smith 
12470971522cSBarry Smith    Level: advanced
12480971522cSBarry Smith 
12490971522cSBarry Smith .seealso: PCFIELDSPLIT
12500971522cSBarry Smith @*/
12517087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12520971522cSBarry Smith {
12534ac538c5SBarry Smith   PetscErrorCode ierr;
12540971522cSBarry Smith 
12550971522cSBarry Smith   PetscFunctionBegin;
12560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
125713e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12584ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12590971522cSBarry Smith   PetscFunctionReturn(0);
12600971522cSBarry Smith }
12610971522cSBarry Smith 
1262e69d4d44SBarry Smith #undef __FUNCT__
1263e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1264e69d4d44SBarry Smith /*@
1265e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1266a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1267e69d4d44SBarry Smith 
1268e69d4d44SBarry Smith     Collective on PC
1269e69d4d44SBarry Smith 
1270e69d4d44SBarry Smith     Input Parameters:
1271e69d4d44SBarry Smith +   pc  - the preconditioner context
1272fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1273084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1274084e4875SJed Brown 
1275e69d4d44SBarry Smith     Options Database:
1276084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1277e69d4d44SBarry Smith 
1278fd1303e9SJungho Lee     Notes:
1279fd1303e9SJungho Lee $    If ptype is
1280fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1281fd1303e9SJungho Lee $             to this function).
1282fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1283fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1284fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1285fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1286fd1303e9SJungho Lee $             preconditioner
1287fd1303e9SJungho Lee 
1288fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1289fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1290fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1291fd1303e9SJungho Lee 
1292fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1293fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1294fd1303e9SJungho Lee 
1295e69d4d44SBarry Smith     Level: intermediate
1296e69d4d44SBarry Smith 
1297fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1298e69d4d44SBarry Smith 
1299e69d4d44SBarry Smith @*/
13007087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1301e69d4d44SBarry Smith {
13024ac538c5SBarry Smith   PetscErrorCode ierr;
1303e69d4d44SBarry Smith 
1304e69d4d44SBarry Smith   PetscFunctionBegin;
13050700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13064ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1307e69d4d44SBarry Smith   PetscFunctionReturn(0);
1308e69d4d44SBarry Smith }
1309e69d4d44SBarry Smith 
1310e69d4d44SBarry Smith EXTERN_C_BEGIN
1311e69d4d44SBarry Smith #undef __FUNCT__
1312e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
13137087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1314e69d4d44SBarry Smith {
1315e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1316084e4875SJed Brown   PetscErrorCode  ierr;
1317e69d4d44SBarry Smith 
1318e69d4d44SBarry Smith   PetscFunctionBegin;
1319084e4875SJed Brown   jac->schurpre = ptype;
1320084e4875SJed Brown   if (pre) {
13216bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1322084e4875SJed Brown     jac->schur_user = pre;
1323084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1324084e4875SJed Brown   }
1325e69d4d44SBarry Smith   PetscFunctionReturn(0);
1326e69d4d44SBarry Smith }
1327e69d4d44SBarry Smith EXTERN_C_END
1328e69d4d44SBarry Smith 
132930ad9308SMatthew Knepley #undef __FUNCT__
1330c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1331ab1df9f5SJed Brown /*@
1332c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1333ab1df9f5SJed Brown 
1334ab1df9f5SJed Brown     Collective on PC
1335ab1df9f5SJed Brown 
1336ab1df9f5SJed Brown     Input Parameters:
1337ab1df9f5SJed Brown +   pc  - the preconditioner context
1338c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1339ab1df9f5SJed Brown 
1340ab1df9f5SJed Brown     Options Database:
1341c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1342ab1df9f5SJed Brown 
1343ab1df9f5SJed Brown 
1344ab1df9f5SJed Brown     Level: intermediate
1345ab1df9f5SJed Brown 
1346ab1df9f5SJed Brown     Notes:
1347ab1df9f5SJed Brown     The FULL factorization is
1348ab1df9f5SJed Brown 
1349ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1350ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1351ab1df9f5SJed Brown 
13526be4592eSBarry 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,
13536be4592eSBarry 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).
1354ab1df9f5SJed Brown 
13556be4592eSBarry 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
13566be4592eSBarry 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
13576be4592eSBarry 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,
13586be4592eSBarry 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.
1359ab1df9f5SJed Brown 
13606be4592eSBarry 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
13616be4592eSBarry 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).
1362ab1df9f5SJed Brown 
1363ab1df9f5SJed Brown     References:
1364ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1365ab1df9f5SJed Brown 
1366ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1367ab1df9f5SJed Brown 
1368ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1369ab1df9f5SJed Brown @*/
1370c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1371ab1df9f5SJed Brown {
1372ab1df9f5SJed Brown   PetscErrorCode ierr;
1373ab1df9f5SJed Brown 
1374ab1df9f5SJed Brown   PetscFunctionBegin;
1375ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1376c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1377ab1df9f5SJed Brown   PetscFunctionReturn(0);
1378ab1df9f5SJed Brown }
1379ab1df9f5SJed Brown 
1380ab1df9f5SJed Brown #undef __FUNCT__
1381c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1382c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1383ab1df9f5SJed Brown {
1384ab1df9f5SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1385ab1df9f5SJed Brown 
1386ab1df9f5SJed Brown   PetscFunctionBegin;
1387ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1388ab1df9f5SJed Brown   PetscFunctionReturn(0);
1389ab1df9f5SJed Brown }
1390ab1df9f5SJed Brown 
1391ab1df9f5SJed Brown #undef __FUNCT__
139230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
139330ad9308SMatthew Knepley /*@C
13948c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
139530ad9308SMatthew Knepley 
139630ad9308SMatthew Knepley    Collective on KSP
139730ad9308SMatthew Knepley 
139830ad9308SMatthew Knepley    Input Parameter:
139930ad9308SMatthew Knepley .  pc - the preconditioner context
140030ad9308SMatthew Knepley 
140130ad9308SMatthew Knepley    Output Parameters:
1402a04f6461SBarry Smith +  A00 - the (0,0) block
1403a04f6461SBarry Smith .  A01 - the (0,1) block
1404a04f6461SBarry Smith .  A10 - the (1,0) block
1405a04f6461SBarry Smith -  A11 - the (1,1) block
140630ad9308SMatthew Knepley 
140730ad9308SMatthew Knepley    Level: advanced
140830ad9308SMatthew Knepley 
140930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
141030ad9308SMatthew Knepley @*/
1411a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
141230ad9308SMatthew Knepley {
141330ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
141430ad9308SMatthew Knepley 
141530ad9308SMatthew Knepley   PetscFunctionBegin;
14160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1417c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1418a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1419a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1420a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1421a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
142230ad9308SMatthew Knepley   PetscFunctionReturn(0);
142330ad9308SMatthew Knepley }
142430ad9308SMatthew Knepley 
142579416396SBarry Smith EXTERN_C_BEGIN
142679416396SBarry Smith #undef __FUNCT__
142779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
14287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
142979416396SBarry Smith {
143079416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1431e69d4d44SBarry Smith   PetscErrorCode ierr;
143279416396SBarry Smith 
143379416396SBarry Smith   PetscFunctionBegin;
143479416396SBarry Smith   jac->type = type;
14353b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
14363b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
14373b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1438e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1439e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1440c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1441e69d4d44SBarry Smith 
14423b224e63SBarry Smith   } else {
14433b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
14443b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1445e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14469e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1447c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
14483b224e63SBarry Smith   }
144979416396SBarry Smith   PetscFunctionReturn(0);
145079416396SBarry Smith }
145179416396SBarry Smith EXTERN_C_END
145279416396SBarry Smith 
145351f519a2SBarry Smith EXTERN_C_BEGIN
145451f519a2SBarry Smith #undef __FUNCT__
145551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
14567087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
145751f519a2SBarry Smith {
145851f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
145951f519a2SBarry Smith 
146051f519a2SBarry Smith   PetscFunctionBegin;
146165e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
146265e19b50SBarry 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);
146351f519a2SBarry Smith   jac->bs = bs;
146451f519a2SBarry Smith   PetscFunctionReturn(0);
146551f519a2SBarry Smith }
146651f519a2SBarry Smith EXTERN_C_END
146751f519a2SBarry Smith 
146879416396SBarry Smith #undef __FUNCT__
146979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1470bc08b0f1SBarry Smith /*@
147179416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
147279416396SBarry Smith 
147379416396SBarry Smith    Collective on PC
147479416396SBarry Smith 
147579416396SBarry Smith    Input Parameter:
147679416396SBarry Smith .  pc - the preconditioner context
147781540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
147879416396SBarry Smith 
147979416396SBarry Smith    Options Database Key:
1480a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
148179416396SBarry Smith 
148279416396SBarry Smith    Level: Developer
148379416396SBarry Smith 
148479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
148579416396SBarry Smith 
148679416396SBarry Smith .seealso: PCCompositeSetType()
148779416396SBarry Smith 
148879416396SBarry Smith @*/
14897087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
149079416396SBarry Smith {
14914ac538c5SBarry Smith   PetscErrorCode ierr;
149279416396SBarry Smith 
149379416396SBarry Smith   PetscFunctionBegin;
14940700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14954ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
149679416396SBarry Smith   PetscFunctionReturn(0);
149779416396SBarry Smith }
149879416396SBarry Smith 
14990971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
15000971522cSBarry Smith /*MC
1501a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1502a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
15030971522cSBarry Smith 
1504edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1505edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
15060971522cSBarry Smith 
1507a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
150869a612a9SBarry Smith          and set the options directly on the resulting KSP object
15090971522cSBarry Smith 
15100971522cSBarry Smith    Level: intermediate
15110971522cSBarry Smith 
151279416396SBarry Smith    Options Database Keys:
151381540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
151481540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
151581540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
151681540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
151781540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1518e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1519435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
152079416396SBarry Smith 
15215d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
15225d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
15235d4c12cdSJungho Lee 
1524c8a0d604SMatthew G Knepley    Notes:
1525c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1526d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1527d32f9abdSBarry Smith 
1528d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1529d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1530d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1531d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1532d32f9abdSBarry Smith 
1533c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1534c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1535c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1536c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1537c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1538a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1539c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1540c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1541c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1542a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1543c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1544c9c6ffaaSJed Brown      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
15455668aaf4SBarry Smith      diag gives
1546c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1547c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
15485668aaf4SBarry 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
1549c8a0d604SMatthew G Knepley $              (  A00   0 )
1550c8a0d604SMatthew G Knepley $              (  A10   S )
1551c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1552c8a0d604SMatthew G Knepley $              ( A00 A01 )
1553c8a0d604SMatthew G Knepley $              (  0   S  )
1554c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1555e69d4d44SBarry Smith 
1556edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1557edf189efSBarry Smith      is used automatically for a second block.
1558edf189efSBarry Smith 
1559ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1560ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1561ff218e97SBarry Smith 
1562ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1563ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1564ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
15650716a85fSBarry Smith 
1566a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
15670971522cSBarry Smith 
15687e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1569e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
15700971522cSBarry Smith M*/
15710971522cSBarry Smith 
15720971522cSBarry Smith EXTERN_C_BEGIN
15730971522cSBarry Smith #undef __FUNCT__
15740971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
15757087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
15760971522cSBarry Smith {
15770971522cSBarry Smith   PetscErrorCode ierr;
15780971522cSBarry Smith   PC_FieldSplit  *jac;
15790971522cSBarry Smith 
15800971522cSBarry Smith   PetscFunctionBegin;
158138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
15820971522cSBarry Smith   jac->bs        = -1;
15830971522cSBarry Smith   jac->nsplits   = 0;
15843e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1585e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1586c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
158751f519a2SBarry Smith 
15880971522cSBarry Smith   pc->data     = (void*)jac;
15890971522cSBarry Smith 
15900971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1591421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
15920971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1593574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
15940971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
15950971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
15960971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
15970971522cSBarry Smith   pc->ops->applyrichardson   = 0;
15980971522cSBarry Smith 
159969a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
160069a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
16010971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
16020971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1603704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1604704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
160579416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
160679416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
160751f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
160851f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
16090971522cSBarry Smith   PetscFunctionReturn(0);
16100971522cSBarry Smith }
16110971522cSBarry Smith EXTERN_C_END
16120971522cSBarry Smith 
16130971522cSBarry Smith 
1614a541d17aSBarry Smith 
1615