xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 939b8a20d39cb19df6ece01ab2377b8403c87da5)
1dba47a55SKris Buschelman 
2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.h" I*/
40971522cSBarry Smith 
5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
7c5d2311dSJed Brown 
80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
90971522cSBarry Smith struct _PC_FieldSplitLink {
1069a612a9SBarry Smith   KSP               ksp;
110971522cSBarry Smith   Vec               x,y;
12db4c96c1SJed Brown   char              *splitname;
130971522cSBarry Smith   PetscInt          nfields;
145d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
151b9fc7fcSBarry Smith   VecScatter        sctx;
165d4c12cdSJungho Lee   IS                is,is_col;
1751f519a2SBarry Smith   PC_FieldSplitLink next,previous;
180971522cSBarry Smith };
190971522cSBarry Smith 
200971522cSBarry Smith typedef struct {
2168dd23aaSBarry Smith   PCCompositeType                    type;
22ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
23ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
24ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
2530ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
2630ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
2779416396SBarry Smith   Vec                                *x,*y,w1,w2;
28519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
29519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3030ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
31ace3abfcSBarry Smith   PetscBool                          issetup;
3230ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3330ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
3430ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
35a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
36084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
37084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
38c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
3930ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4097bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4163ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
42c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
430971522cSBarry Smith } PC_FieldSplit;
440971522cSBarry Smith 
4516913363SBarry Smith /*
4616913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
4716913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
4816913363SBarry Smith    PC you could change this.
4916913363SBarry Smith */
50084e4875SJed Brown 
51e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
52084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
53084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
54084e4875SJed Brown {
55084e4875SJed Brown   switch (jac->schurpre) {
56084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
57084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
58084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
59084e4875SJed Brown     default:
60084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
61084e4875SJed Brown   }
62084e4875SJed Brown }
63084e4875SJed Brown 
64084e4875SJed Brown 
650971522cSBarry Smith #undef __FUNCT__
660971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
670971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
680971522cSBarry Smith {
690971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
700971522cSBarry Smith   PetscErrorCode    ierr;
71ace3abfcSBarry Smith   PetscBool         iascii;
720971522cSBarry Smith   PetscInt          i,j;
735a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
740971522cSBarry Smith 
750971522cSBarry Smith   PetscFunctionBegin;
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;
225e7c4fc90SDmitry Karpeev       DM          *dms;
226e7c4fc90SDmitry Karpeev       ierr = DMCreateDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);    CHKERRQ(ierr);
2277b62db95SJungho Lee       for(f = 0; f < numFields; ++f) {
2287b62db95SJungho Lee         ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
2297b62db95SJungho Lee         ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
2307b62db95SJungho Lee         ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
2317b62db95SJungho Lee       }
2327b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
2337b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
234e7c4fc90SDmitry Karpeev       if(dms) {
2358b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2367b62db95SJungho Lee         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2377b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
2387b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
239e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]); CHKERRQ(ierr);
2402fa5ba8aSJed Brown         }
2417b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
2428b8307b2SJed Brown       }
24366ffff09SJed Brown     } else {
244521d7252SBarry Smith       if (jac->bs <= 0) {
245704ba839SBarry Smith         if (pc->pmat) {
246521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
247704ba839SBarry Smith         } else {
248704ba839SBarry Smith           jac->bs = 1;
249704ba839SBarry Smith         }
250521d7252SBarry Smith       }
251d32f9abdSBarry Smith 
252acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2536ce1633cSBarry Smith       if (stokes) {
2546ce1633cSBarry Smith         IS       zerodiags,rest;
2556ce1633cSBarry Smith         PetscInt nmin,nmax;
2566ce1633cSBarry Smith 
2576ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2586ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2596ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2606ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2616ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
262fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
263fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2646ce1633cSBarry Smith       } else {
265d32f9abdSBarry Smith         if (!flg) {
266d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
267d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2686c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2696c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
270d32f9abdSBarry Smith         }
2716c924f48SJed Brown         if (flg || !jac->splitdefined) {
272d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
273db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2746c924f48SJed Brown             char splitname[8];
2756c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2765d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
27779416396SBarry Smith           }
2785d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
279521d7252SBarry Smith         }
28066ffff09SJed Brown       }
2816ce1633cSBarry Smith     }
282edf189efSBarry Smith   } else if (jac->nsplits == 1) {
283edf189efSBarry Smith     if (ilink->is) {
284edf189efSBarry Smith       IS       is2;
285edf189efSBarry Smith       PetscInt nmin,nmax;
286edf189efSBarry Smith 
287edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
288edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
289db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
290fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
2917b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
29263ec74ffSBarry Smith   } else if (jac->reset) {
29363ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
294d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
295d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
296d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
297d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
298d0af7cd3SBarry Smith       PetscBool dmcomposite;
299d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
300d0af7cd3SBarry Smith       if (dmcomposite) {
301d0af7cd3SBarry Smith         PetscInt nDM;
302d0af7cd3SBarry Smith         IS       *fields;
3037b62db95SJungho Lee         DM       *dms;
304d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
305d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
306d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3077b62db95SJungho Lee         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3087b62db95SJungho Lee         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
309d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3107b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3117b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
312d0af7cd3SBarry Smith           ilink->is = fields[i];
313d0af7cd3SBarry Smith           ilink     = ilink->next;
314edf189efSBarry Smith         }
315d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3167b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
317d0af7cd3SBarry Smith       }
318d0af7cd3SBarry Smith     } else {
319d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
320d0af7cd3SBarry Smith       if (stokes) {
321d0af7cd3SBarry Smith         IS       zerodiags,rest;
322d0af7cd3SBarry Smith         PetscInt nmin,nmax;
323d0af7cd3SBarry Smith 
324d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
325d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
326d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
32720b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
32820b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
329d0af7cd3SBarry Smith         ilink->is       = rest;
330d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
33120b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
332d0af7cd3SBarry Smith     }
333d0af7cd3SBarry Smith   }
334d0af7cd3SBarry Smith 
3357b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
33669a612a9SBarry Smith   PetscFunctionReturn(0);
33769a612a9SBarry Smith }
33869a612a9SBarry Smith 
33969a612a9SBarry Smith #undef __FUNCT__
34069a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
34169a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
34269a612a9SBarry Smith {
34369a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
34469a612a9SBarry Smith   PetscErrorCode    ierr;
3455a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
346cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
34769a612a9SBarry Smith   MatStructure      flag = pc->flag;
3485f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
34969a612a9SBarry Smith 
35069a612a9SBarry Smith   PetscFunctionBegin;
35169a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
35297bbdb24SBarry Smith   nsplit = jac->nsplits;
3535a9f2f41SSatish Balay   ilink  = jac->head;
35497bbdb24SBarry Smith 
35597bbdb24SBarry Smith   /* get the matrices for each split */
356704ba839SBarry Smith   if (!jac->issetup) {
3571b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
35897bbdb24SBarry Smith 
359704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
360704ba839SBarry Smith 
3615d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
36251f519a2SBarry Smith     bs     = jac->bs;
36397bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
364cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3651b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3661b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3671b9fc7fcSBarry Smith       if (jac->defaultsplit) {
368704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
3695f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
370704ba839SBarry Smith       } else if (!ilink->is) {
371ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3725f4ab4e1SJungho Lee           PetscInt   *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
3735a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3745f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
3751b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3761b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3771b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
3785f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
37997bbdb24SBarry Smith             }
38097bbdb24SBarry Smith           }
381d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
3825f4ab4e1SJungho Lee           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
383ccb205f8SBarry Smith         } else {
384704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
3855f4ab4e1SJungho Lee           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
386ccb205f8SBarry Smith        }
3873e197d65SBarry Smith       }
388704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
3895f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
3905f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
391704ba839SBarry Smith       ilink = ilink->next;
3921b9fc7fcSBarry Smith     }
3931b9fc7fcSBarry Smith   }
3941b9fc7fcSBarry Smith 
395704ba839SBarry Smith   ilink  = jac->head;
39697bbdb24SBarry Smith   if (!jac->pmat) {
397cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
398cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3995f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
400704ba839SBarry Smith       ilink = ilink->next;
401cf502942SBarry Smith     }
40297bbdb24SBarry Smith   } else {
403cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4045f4ab4e1SJungho Lee       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
405704ba839SBarry Smith       ilink = ilink->next;
406cf502942SBarry Smith     }
40797bbdb24SBarry Smith   }
408519d70e2SJed Brown   if (jac->realdiagonal) {
409519d70e2SJed Brown     ilink = jac->head;
410519d70e2SJed Brown     if (!jac->mat) {
411519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
412519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4135f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
414519d70e2SJed Brown         ilink = ilink->next;
415519d70e2SJed Brown       }
416519d70e2SJed Brown     } else {
417519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
4185f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
419519d70e2SJed Brown         ilink = ilink->next;
420519d70e2SJed Brown       }
421519d70e2SJed Brown     }
422519d70e2SJed Brown   } else {
423519d70e2SJed Brown     jac->mat = jac->pmat;
424519d70e2SJed Brown   }
42597bbdb24SBarry Smith 
4266c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
42768dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
42868dd23aaSBarry Smith     ilink  = jac->head;
42968dd23aaSBarry Smith     if (!jac->Afield) {
43068dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
43168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4324aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
43368dd23aaSBarry Smith         ilink = ilink->next;
43468dd23aaSBarry Smith       }
43568dd23aaSBarry Smith     } else {
43668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4374aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
43868dd23aaSBarry Smith         ilink = ilink->next;
43968dd23aaSBarry Smith       }
44068dd23aaSBarry Smith     }
44168dd23aaSBarry Smith   }
44268dd23aaSBarry Smith 
4433b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4443b224e63SBarry Smith     IS       ccis;
4454aa3045dSJed Brown     PetscInt rstart,rend;
446093c86ffSJed Brown     char     lscname[256];
447093c86ffSJed Brown     PetscObject LSC_L;
448e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
44968dd23aaSBarry Smith 
450e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
451e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
452e6cab6aaSJed Brown 
4533b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4543b224e63SBarry Smith     if (jac->schur) {
4553b224e63SBarry Smith       ilink = jac->head;
45649bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4574aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
458fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4593b224e63SBarry Smith       ilink = ilink->next;
46049bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4614aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
462fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
463519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
464084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4653b224e63SBarry Smith 
4663b224e63SBarry Smith      } else {
4671cee3971SBarry Smith       KSP ksp;
4686c924f48SJed Brown       char schurprefix[256];
4693b224e63SBarry Smith 
470a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4713b224e63SBarry Smith       ilink = jac->head;
47249bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4734aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
474fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4753b224e63SBarry Smith       ilink = ilink->next;
47649bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
4774aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
478fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4797d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
480519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
481a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4821cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
483e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
484a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
485a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
48620b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
48720b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4887b62db95SJungho Lee       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
4893b224e63SBarry Smith 
4903b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4919005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4921cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
493084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
494084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
495084e4875SJed Brown         PC pc;
496cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
497084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
498084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
499e69d4d44SBarry Smith       }
5006c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5016c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5023b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
50320b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
50420b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5053b224e63SBarry Smith 
5063b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5073b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5083b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5093b224e63SBarry Smith       ilink = jac->head;
5103b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5113b224e63SBarry Smith       ilink = ilink->next;
5123b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5133b224e63SBarry Smith     }
514093c86ffSJed Brown 
515093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
516093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
517093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
518093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
519093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
520093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
521093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
522093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
523093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
5243b224e63SBarry Smith   } else {
52597bbdb24SBarry Smith     /* set up the individual PCs */
52697bbdb24SBarry Smith     i    = 0;
5275a9f2f41SSatish Balay     ilink = jac->head;
5285a9f2f41SSatish Balay     while (ilink) {
529519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5303b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
531c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5325a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
53397bbdb24SBarry Smith       i++;
5345a9f2f41SSatish Balay       ilink = ilink->next;
5350971522cSBarry Smith     }
53697bbdb24SBarry Smith 
53797bbdb24SBarry Smith     /* create work vectors for each split */
5381b9fc7fcSBarry Smith     if (!jac->x) {
53997bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5405a9f2f41SSatish Balay       ilink = jac->head;
54197bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
542906ed7ccSBarry Smith         Vec *vl,*vr;
543906ed7ccSBarry Smith 
544906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
545906ed7ccSBarry Smith         ilink->x  = *vr;
546906ed7ccSBarry Smith         ilink->y  = *vl;
547906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
548906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5495a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5505a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5515a9f2f41SSatish Balay         ilink     = ilink->next;
55297bbdb24SBarry Smith       }
5533b224e63SBarry Smith     }
5543b224e63SBarry Smith   }
5553b224e63SBarry Smith 
5563b224e63SBarry Smith 
557d0f46423SBarry Smith   if (!jac->head->sctx) {
5583b224e63SBarry Smith     Vec xtmp;
5593b224e63SBarry Smith 
56079416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5611b9fc7fcSBarry Smith 
5625a9f2f41SSatish Balay     ilink = jac->head;
5631b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5641b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
565704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5665a9f2f41SSatish Balay       ilink = ilink->next;
5671b9fc7fcSBarry Smith     }
5686bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5691b9fc7fcSBarry Smith   }
570c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5710971522cSBarry Smith   PetscFunctionReturn(0);
5720971522cSBarry Smith }
5730971522cSBarry Smith 
5745a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
575ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
576ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5775a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
578ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
579ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
58079416396SBarry Smith 
5810971522cSBarry Smith #undef __FUNCT__
5823b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5833b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5843b224e63SBarry Smith {
5853b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5863b224e63SBarry Smith   PetscErrorCode    ierr;
5873b224e63SBarry Smith   KSP               ksp;
5883b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5893b224e63SBarry Smith 
5903b224e63SBarry Smith   PetscFunctionBegin;
5913b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5923b224e63SBarry Smith 
593c5d2311dSJed Brown   switch (jac->schurfactorization) {
594c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
595a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
596c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
597c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
599c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
600c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
601c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
602c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
603c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
604c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
605c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
606c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
607c5d2311dSJed Brown       break;
608c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
609a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
610c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
611c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
615c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
616c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
617c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
618c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
619c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
620c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
621c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
622c5d2311dSJed Brown       break;
623c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
624a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
625c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
628c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
629c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
630c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
631c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
632c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
633c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
634c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
635c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
636c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
637c5d2311dSJed Brown       break;
638c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_FULL:
639a04f6461SBarry 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 */
6403b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6413b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6423b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6433b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6443b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6453b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6463b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6473b224e63SBarry Smith 
6483b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6493b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6503b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6513b224e63SBarry Smith 
6523b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6533b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6543b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6553b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6563b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
657c5d2311dSJed Brown   }
6583b224e63SBarry Smith   PetscFunctionReturn(0);
6593b224e63SBarry Smith }
6603b224e63SBarry Smith 
6613b224e63SBarry Smith #undef __FUNCT__
6620971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6630971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6640971522cSBarry Smith {
6650971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6660971522cSBarry Smith   PetscErrorCode    ierr;
6675a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
668*939b8a20SBarry Smith   PetscInt          cnt,bs;
6690971522cSBarry Smith 
6700971522cSBarry Smith   PetscFunctionBegin;
67151f519a2SBarry Smith   CHKMEMQ;
67279416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6731b9fc7fcSBarry Smith     if (jac->defaultsplit) {
674*939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
675*939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
676*939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
677*939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
6780971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6795a9f2f41SSatish Balay       while (ilink) {
6805a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6815a9f2f41SSatish Balay         ilink = ilink->next;
6820971522cSBarry Smith       }
6830971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6841b9fc7fcSBarry Smith     } else {
685efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6865a9f2f41SSatish Balay       while (ilink) {
6875a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6885a9f2f41SSatish Balay         ilink = ilink->next;
6891b9fc7fcSBarry Smith       }
6901b9fc7fcSBarry Smith     }
69116913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
69279416396SBarry Smith     if (!jac->w1) {
69379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
69479416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
69579416396SBarry Smith     }
696efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6975a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6983e197d65SBarry Smith     cnt = 1;
6995a9f2f41SSatish Balay     while (ilink->next) {
7005a9f2f41SSatish Balay       ilink = ilink->next;
7013e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7023e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7033e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7043e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7053e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7063e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7073e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7083e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7093e197d65SBarry Smith     }
71051f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
71111755939SBarry Smith       cnt -= 2;
71251f519a2SBarry Smith       while (ilink->previous) {
71351f519a2SBarry Smith         ilink = ilink->previous;
71411755939SBarry Smith         /* compute the residual only over the part of the vector needed */
71511755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
71611755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
71711755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71811755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71911755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
72011755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72111755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72251f519a2SBarry Smith       }
72311755939SBarry Smith     }
72465e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
72551f519a2SBarry Smith   CHKMEMQ;
7260971522cSBarry Smith   PetscFunctionReturn(0);
7270971522cSBarry Smith }
7280971522cSBarry Smith 
729421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
730ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
731ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
732421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
733ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
734ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
735421e10b8SBarry Smith 
736421e10b8SBarry Smith #undef __FUNCT__
7378c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
738421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
739421e10b8SBarry Smith {
740421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
741421e10b8SBarry Smith   PetscErrorCode    ierr;
742421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
743*939b8a20SBarry Smith   PetscInt          bs;
744421e10b8SBarry Smith 
745421e10b8SBarry Smith   PetscFunctionBegin;
746421e10b8SBarry Smith   CHKMEMQ;
747421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
748421e10b8SBarry Smith     if (jac->defaultsplit) {
749*939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
750*939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
751*939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
752*939b8a20SBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
753421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
754421e10b8SBarry Smith       while (ilink) {
755421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
756421e10b8SBarry Smith 	ilink = ilink->next;
757421e10b8SBarry Smith       }
758421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
759421e10b8SBarry Smith     } else {
760421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
761421e10b8SBarry Smith       while (ilink) {
762421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
763421e10b8SBarry Smith 	ilink = ilink->next;
764421e10b8SBarry Smith       }
765421e10b8SBarry Smith     }
766421e10b8SBarry Smith   } else {
767421e10b8SBarry Smith     if (!jac->w1) {
768421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
769421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
770421e10b8SBarry Smith     }
771421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
772421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
773421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
774421e10b8SBarry Smith       while (ilink->next) {
775421e10b8SBarry Smith         ilink = ilink->next;
7769989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
777421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
778421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
779421e10b8SBarry Smith       }
780421e10b8SBarry Smith       while (ilink->previous) {
781421e10b8SBarry Smith         ilink = ilink->previous;
7829989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
783421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
784421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
785421e10b8SBarry Smith       }
786421e10b8SBarry Smith     } else {
787421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
788421e10b8SBarry Smith 	ilink = ilink->next;
789421e10b8SBarry Smith       }
790421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
791421e10b8SBarry Smith       while (ilink->previous) {
792421e10b8SBarry Smith 	ilink = ilink->previous;
7939989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
794421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
795421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
796421e10b8SBarry Smith       }
797421e10b8SBarry Smith     }
798421e10b8SBarry Smith   }
799421e10b8SBarry Smith   CHKMEMQ;
800421e10b8SBarry Smith   PetscFunctionReturn(0);
801421e10b8SBarry Smith }
802421e10b8SBarry Smith 
8030971522cSBarry Smith #undef __FUNCT__
804574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
805574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8060971522cSBarry Smith {
8070971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8080971522cSBarry Smith   PetscErrorCode    ierr;
8095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8100971522cSBarry Smith 
8110971522cSBarry Smith   PetscFunctionBegin;
8125a9f2f41SSatish Balay   while (ilink) {
813574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
814fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
815fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
816fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
817fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
818b5787286SJed Brown     ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
8195a9f2f41SSatish Balay     next = ilink->next;
8205a9f2f41SSatish Balay     ilink = next;
8210971522cSBarry Smith   }
82205b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
823574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
824574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
825574deadeSBarry Smith   } else if (jac->mat) {
826574deadeSBarry Smith     jac->mat = PETSC_NULL;
827574deadeSBarry Smith   }
82897bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
82968dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8306bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8316bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8326bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8336bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8346bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8356bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8366bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
83763ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
838574deadeSBarry Smith   PetscFunctionReturn(0);
839574deadeSBarry Smith }
840574deadeSBarry Smith 
841574deadeSBarry Smith #undef __FUNCT__
842574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
843574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
844574deadeSBarry Smith {
845574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
846574deadeSBarry Smith   PetscErrorCode    ierr;
847574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
848574deadeSBarry Smith 
849574deadeSBarry Smith   PetscFunctionBegin;
850574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
851574deadeSBarry Smith   while (ilink) {
8526bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
853574deadeSBarry Smith     next = ilink->next;
854574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
855574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
8565d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
857574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
858574deadeSBarry Smith     ilink = next;
859574deadeSBarry Smith   }
860574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
861c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8627b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
8637b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
8647b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
8657b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
8667b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
867ab1df9f5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr);
868c9c6ffaaSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr);
8690971522cSBarry Smith   PetscFunctionReturn(0);
8700971522cSBarry Smith }
8710971522cSBarry Smith 
8720971522cSBarry Smith #undef __FUNCT__
8730971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8740971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8750971522cSBarry Smith {
8761b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8776c924f48SJed Brown   PetscInt        bs;
878bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8799dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8803b224e63SBarry Smith   PCCompositeType ctype;
881e7c4fc90SDmitry Karpeev   DM              ddm;
882e7c4fc90SDmitry Karpeev   char            ddm_name[1025];
8831b9fc7fcSBarry Smith 
8840971522cSBarry Smith   PetscFunctionBegin;
8851b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
886e7c4fc90SDmitry Karpeev   if(pc->dm) {
887e7c4fc90SDmitry Karpeev     /* Allow the user to request a decomposition DM by name */
888e7c4fc90SDmitry Karpeev     ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
889e7c4fc90SDmitry Karpeev     ierr = PetscOptionsString("-pc_fieldsplit_decomposition_dm", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
890e7c4fc90SDmitry Karpeev     if(flg) {
891e7c4fc90SDmitry Karpeev       ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
892e7c4fc90SDmitry Karpeev       if(!ddm) {
893e7c4fc90SDmitry Karpeev         SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
894e7c4fc90SDmitry Karpeev       }
895e7c4fc90SDmitry Karpeev       ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
896e7c4fc90SDmitry Karpeev       ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
897e7c4fc90SDmitry Karpeev     }
898e7c4fc90SDmitry Karpeev   }
899acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
90051f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
90151f519a2SBarry Smith   if (flg) {
90251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
90351f519a2SBarry Smith   }
904704ba839SBarry Smith 
905435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
906c0adfefeSBarry Smith   if (stokes) {
907c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
908c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
909c0adfefeSBarry Smith   }
910c0adfefeSBarry Smith 
9113b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
9123b224e63SBarry Smith   if (flg) {
9133b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
9143b224e63SBarry Smith   }
915c30613efSMatthew Knepley   /* Only setup fields once */
916c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
917d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
918d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
9196c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9206c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
921d32f9abdSBarry Smith   }
922c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
923c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
924c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
925c9c6ffaaSJed 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);
926084e4875SJed 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);
927c5d2311dSJed Brown   }
9281b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9290971522cSBarry Smith   PetscFunctionReturn(0);
9300971522cSBarry Smith }
9310971522cSBarry Smith 
9320971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9330971522cSBarry Smith 
9340971522cSBarry Smith EXTERN_C_BEGIN
9350971522cSBarry Smith #undef __FUNCT__
9360971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9375d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
9380971522cSBarry Smith {
93997bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9400971522cSBarry Smith   PetscErrorCode    ierr;
9415a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
94269a612a9SBarry Smith   char              prefix[128];
9435d4c12cdSJungho Lee   PetscInt          i;
9440971522cSBarry Smith 
9450971522cSBarry Smith   PetscFunctionBegin;
9466c924f48SJed Brown   if (jac->splitdefined) {
9476c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9486c924f48SJed Brown     PetscFunctionReturn(0);
9496c924f48SJed Brown   }
95051f519a2SBarry Smith   for (i=0; i<n; i++) {
951e32f2f54SBarry 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);
952e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
95351f519a2SBarry Smith   }
954704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
955a04f6461SBarry Smith   if (splitname) {
956db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
957a04f6461SBarry Smith   } else {
958a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
959a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
960a04f6461SBarry Smith   }
961704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9625a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9635d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
9645d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
9655a9f2f41SSatish Balay   ilink->nfields = n;
9665a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9677adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9681cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9695a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9709005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
97169a612a9SBarry Smith 
972a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9735a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9740971522cSBarry Smith 
9750971522cSBarry Smith   if (!next) {
9765a9f2f41SSatish Balay     jac->head       = ilink;
97751f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9780971522cSBarry Smith   } else {
9790971522cSBarry Smith     while (next->next) {
9800971522cSBarry Smith       next = next->next;
9810971522cSBarry Smith     }
9825a9f2f41SSatish Balay     next->next      = ilink;
98351f519a2SBarry Smith     ilink->previous = next;
9840971522cSBarry Smith   }
9850971522cSBarry Smith   jac->nsplits++;
9860971522cSBarry Smith   PetscFunctionReturn(0);
9870971522cSBarry Smith }
9880971522cSBarry Smith EXTERN_C_END
9890971522cSBarry Smith 
990e69d4d44SBarry Smith EXTERN_C_BEGIN
991e69d4d44SBarry Smith #undef __FUNCT__
992e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9937087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
994e69d4d44SBarry Smith {
995e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
996e69d4d44SBarry Smith   PetscErrorCode ierr;
997e69d4d44SBarry Smith 
998e69d4d44SBarry Smith   PetscFunctionBegin;
999519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1000e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1001e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
100213e0d083SBarry Smith   if (n) *n = jac->nsplits;
1003e69d4d44SBarry Smith   PetscFunctionReturn(0);
1004e69d4d44SBarry Smith }
1005e69d4d44SBarry Smith EXTERN_C_END
10060971522cSBarry Smith 
10070971522cSBarry Smith EXTERN_C_BEGIN
10080971522cSBarry Smith #undef __FUNCT__
100969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
10107087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
10110971522cSBarry Smith {
10120971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10130971522cSBarry Smith   PetscErrorCode    ierr;
10140971522cSBarry Smith   PetscInt          cnt = 0;
10155a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
10160971522cSBarry Smith 
10170971522cSBarry Smith   PetscFunctionBegin;
10185d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
10195a9f2f41SSatish Balay   while (ilink) {
10205a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
10215a9f2f41SSatish Balay     ilink = ilink->next;
10220971522cSBarry Smith   }
10235d480477SMatthew 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);
102413e0d083SBarry Smith   if (n) *n = jac->nsplits;
10250971522cSBarry Smith   PetscFunctionReturn(0);
10260971522cSBarry Smith }
10270971522cSBarry Smith EXTERN_C_END
10280971522cSBarry Smith 
1029704ba839SBarry Smith EXTERN_C_BEGIN
1030704ba839SBarry Smith #undef __FUNCT__
1031704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10327087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1033704ba839SBarry Smith {
1034704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1035704ba839SBarry Smith   PetscErrorCode    ierr;
1036704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1037704ba839SBarry Smith   char              prefix[128];
1038704ba839SBarry Smith 
1039704ba839SBarry Smith   PetscFunctionBegin;
10406c924f48SJed Brown   if (jac->splitdefined) {
10416c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10426c924f48SJed Brown     PetscFunctionReturn(0);
10436c924f48SJed Brown   }
104416913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1045a04f6461SBarry Smith   if (splitname) {
1046db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1047a04f6461SBarry Smith   } else {
1048b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1049b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1050a04f6461SBarry Smith   }
10511ab39975SBarry Smith   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1052b5787286SJed Brown   ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
1053b5787286SJed Brown   ilink->is      = is;
1054b5787286SJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1055b5787286SJed Brown   ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1056b5787286SJed Brown   ilink->is_col  = is;
1057704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1058704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10591cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1060704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10619005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1062704ba839SBarry Smith 
1063a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1064704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1065704ba839SBarry Smith 
1066704ba839SBarry Smith   if (!next) {
1067704ba839SBarry Smith     jac->head       = ilink;
1068704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1069704ba839SBarry Smith   } else {
1070704ba839SBarry Smith     while (next->next) {
1071704ba839SBarry Smith       next = next->next;
1072704ba839SBarry Smith     }
1073704ba839SBarry Smith     next->next      = ilink;
1074704ba839SBarry Smith     ilink->previous = next;
1075704ba839SBarry Smith   }
1076704ba839SBarry Smith   jac->nsplits++;
1077704ba839SBarry Smith 
1078704ba839SBarry Smith   PetscFunctionReturn(0);
1079704ba839SBarry Smith }
1080704ba839SBarry Smith EXTERN_C_END
1081704ba839SBarry Smith 
10820971522cSBarry Smith #undef __FUNCT__
10830971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10840971522cSBarry Smith /*@
10850971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10860971522cSBarry Smith 
1087ad4df100SBarry Smith     Logically Collective on PC
10880971522cSBarry Smith 
10890971522cSBarry Smith     Input Parameters:
10900971522cSBarry Smith +   pc  - the preconditioner context
1091a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10920971522cSBarry Smith .   n - the number of fields in this split
1093db4c96c1SJed Brown -   fields - the fields in this split
10940971522cSBarry Smith 
10950971522cSBarry Smith     Level: intermediate
10960971522cSBarry Smith 
1097d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1098d32f9abdSBarry Smith 
1099d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1100d32f9abdSBarry 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
1101d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1102d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1103d32f9abdSBarry Smith 
1104db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1105db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1106db4c96c1SJed Brown 
11075d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
11085d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
11095d4c12cdSJungho Lee      available when this routine is called.
11105d4c12cdSJungho Lee 
1111d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
11120971522cSBarry Smith 
11130971522cSBarry Smith @*/
11145d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
11150971522cSBarry Smith {
11164ac538c5SBarry Smith   PetscErrorCode ierr;
11170971522cSBarry Smith 
11180971522cSBarry Smith   PetscFunctionBegin;
11190700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1120db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1121db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1122db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
11235d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
11240971522cSBarry Smith   PetscFunctionReturn(0);
11250971522cSBarry Smith }
11260971522cSBarry Smith 
11270971522cSBarry Smith #undef __FUNCT__
1128704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1129704ba839SBarry Smith /*@
1130704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1131704ba839SBarry Smith 
1132ad4df100SBarry Smith     Logically Collective on PC
1133704ba839SBarry Smith 
1134704ba839SBarry Smith     Input Parameters:
1135704ba839SBarry Smith +   pc  - the preconditioner context
1136a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1137db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1138704ba839SBarry Smith 
1139d32f9abdSBarry Smith 
1140a6ffb8dbSJed Brown     Notes:
1141a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1142a6ffb8dbSJed Brown 
1143db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1144db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1145d32f9abdSBarry Smith 
1146704ba839SBarry Smith     Level: intermediate
1147704ba839SBarry Smith 
1148704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1149704ba839SBarry Smith 
1150704ba839SBarry Smith @*/
11517087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1152704ba839SBarry Smith {
11534ac538c5SBarry Smith   PetscErrorCode ierr;
1154704ba839SBarry Smith 
1155704ba839SBarry Smith   PetscFunctionBegin;
11560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11577b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1158db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11594ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1160704ba839SBarry Smith   PetscFunctionReturn(0);
1161704ba839SBarry Smith }
1162704ba839SBarry Smith 
1163704ba839SBarry Smith #undef __FUNCT__
116457a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
116557a9adfeSMatthew G Knepley /*@
116657a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
116757a9adfeSMatthew G Knepley 
116857a9adfeSMatthew G Knepley     Logically Collective on PC
116957a9adfeSMatthew G Knepley 
117057a9adfeSMatthew G Knepley     Input Parameters:
117157a9adfeSMatthew G Knepley +   pc  - the preconditioner context
117257a9adfeSMatthew G Knepley -   splitname - name of this split
117357a9adfeSMatthew G Knepley 
117457a9adfeSMatthew G Knepley     Output Parameter:
117557a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
117657a9adfeSMatthew G Knepley 
117757a9adfeSMatthew G Knepley     Level: intermediate
117857a9adfeSMatthew G Knepley 
117957a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
118057a9adfeSMatthew G Knepley 
118157a9adfeSMatthew G Knepley @*/
118257a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
118357a9adfeSMatthew G Knepley {
118457a9adfeSMatthew G Knepley   PetscErrorCode ierr;
118557a9adfeSMatthew G Knepley 
118657a9adfeSMatthew G Knepley   PetscFunctionBegin;
118757a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
118857a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
118957a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
119057a9adfeSMatthew G Knepley   {
119157a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
119257a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
119357a9adfeSMatthew G Knepley     PetscBool         found;
119457a9adfeSMatthew G Knepley 
119557a9adfeSMatthew G Knepley     *is = PETSC_NULL;
119657a9adfeSMatthew G Knepley     while(ilink) {
119757a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
119857a9adfeSMatthew G Knepley       if (found) {
119957a9adfeSMatthew G Knepley         *is = ilink->is;
120057a9adfeSMatthew G Knepley         break;
120157a9adfeSMatthew G Knepley       }
120257a9adfeSMatthew G Knepley       ilink = ilink->next;
120357a9adfeSMatthew G Knepley     }
120457a9adfeSMatthew G Knepley   }
120557a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
120657a9adfeSMatthew G Knepley }
120757a9adfeSMatthew G Knepley 
120857a9adfeSMatthew G Knepley #undef __FUNCT__
120951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
121051f519a2SBarry Smith /*@
121151f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
121251f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
121351f519a2SBarry Smith 
1214ad4df100SBarry Smith     Logically Collective on PC
121551f519a2SBarry Smith 
121651f519a2SBarry Smith     Input Parameters:
121751f519a2SBarry Smith +   pc  - the preconditioner context
121851f519a2SBarry Smith -   bs - the block size
121951f519a2SBarry Smith 
122051f519a2SBarry Smith     Level: intermediate
122151f519a2SBarry Smith 
122251f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
122351f519a2SBarry Smith 
122451f519a2SBarry Smith @*/
12257087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
122651f519a2SBarry Smith {
12274ac538c5SBarry Smith   PetscErrorCode ierr;
122851f519a2SBarry Smith 
122951f519a2SBarry Smith   PetscFunctionBegin;
12300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1231c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
12324ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
123351f519a2SBarry Smith   PetscFunctionReturn(0);
123451f519a2SBarry Smith }
123551f519a2SBarry Smith 
123651f519a2SBarry Smith #undef __FUNCT__
123769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12380971522cSBarry Smith /*@C
123969a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12400971522cSBarry Smith 
124169a612a9SBarry Smith    Collective on KSP
12420971522cSBarry Smith 
12430971522cSBarry Smith    Input Parameter:
12440971522cSBarry Smith .  pc - the preconditioner context
12450971522cSBarry Smith 
12460971522cSBarry Smith    Output Parameters:
124713e0d083SBarry Smith +  n - the number of splits
124869a612a9SBarry Smith -  pc - the array of KSP contexts
12490971522cSBarry Smith 
12500971522cSBarry Smith    Note:
1251d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1252d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12530971522cSBarry Smith 
125469a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12550971522cSBarry Smith 
12560971522cSBarry Smith    Level: advanced
12570971522cSBarry Smith 
12580971522cSBarry Smith .seealso: PCFIELDSPLIT
12590971522cSBarry Smith @*/
12607087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12610971522cSBarry Smith {
12624ac538c5SBarry Smith   PetscErrorCode ierr;
12630971522cSBarry Smith 
12640971522cSBarry Smith   PetscFunctionBegin;
12650700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
126613e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12674ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12680971522cSBarry Smith   PetscFunctionReturn(0);
12690971522cSBarry Smith }
12700971522cSBarry Smith 
1271e69d4d44SBarry Smith #undef __FUNCT__
1272e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1273e69d4d44SBarry Smith /*@
1274e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1275a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1276e69d4d44SBarry Smith 
1277e69d4d44SBarry Smith     Collective on PC
1278e69d4d44SBarry Smith 
1279e69d4d44SBarry Smith     Input Parameters:
1280e69d4d44SBarry Smith +   pc  - the preconditioner context
1281fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1282084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1283084e4875SJed Brown 
1284e69d4d44SBarry Smith     Options Database:
1285084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1286e69d4d44SBarry Smith 
1287fd1303e9SJungho Lee     Notes:
1288fd1303e9SJungho Lee $    If ptype is
1289fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1290fd1303e9SJungho Lee $             to this function).
1291fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1292fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1293fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1294fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1295fd1303e9SJungho Lee $             preconditioner
1296fd1303e9SJungho Lee 
1297fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1298fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1299fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1300fd1303e9SJungho Lee 
1301fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1302fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1303fd1303e9SJungho Lee 
1304e69d4d44SBarry Smith     Level: intermediate
1305e69d4d44SBarry Smith 
1306fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1307e69d4d44SBarry Smith 
1308e69d4d44SBarry Smith @*/
13097087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1310e69d4d44SBarry Smith {
13114ac538c5SBarry Smith   PetscErrorCode ierr;
1312e69d4d44SBarry Smith 
1313e69d4d44SBarry Smith   PetscFunctionBegin;
13140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13154ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1316e69d4d44SBarry Smith   PetscFunctionReturn(0);
1317e69d4d44SBarry Smith }
1318e69d4d44SBarry Smith 
1319e69d4d44SBarry Smith EXTERN_C_BEGIN
1320e69d4d44SBarry Smith #undef __FUNCT__
1321e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
13227087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1323e69d4d44SBarry Smith {
1324e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1325084e4875SJed Brown   PetscErrorCode  ierr;
1326e69d4d44SBarry Smith 
1327e69d4d44SBarry Smith   PetscFunctionBegin;
1328084e4875SJed Brown   jac->schurpre = ptype;
1329084e4875SJed Brown   if (pre) {
13306bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1331084e4875SJed Brown     jac->schur_user = pre;
1332084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1333084e4875SJed Brown   }
1334e69d4d44SBarry Smith   PetscFunctionReturn(0);
1335e69d4d44SBarry Smith }
1336e69d4d44SBarry Smith EXTERN_C_END
1337e69d4d44SBarry Smith 
133830ad9308SMatthew Knepley #undef __FUNCT__
1339c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1340ab1df9f5SJed Brown /*@
1341c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1342ab1df9f5SJed Brown 
1343ab1df9f5SJed Brown     Collective on PC
1344ab1df9f5SJed Brown 
1345ab1df9f5SJed Brown     Input Parameters:
1346ab1df9f5SJed Brown +   pc  - the preconditioner context
1347c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1348ab1df9f5SJed Brown 
1349ab1df9f5SJed Brown     Options Database:
1350c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1351ab1df9f5SJed Brown 
1352ab1df9f5SJed Brown 
1353ab1df9f5SJed Brown     Level: intermediate
1354ab1df9f5SJed Brown 
1355ab1df9f5SJed Brown     Notes:
1356ab1df9f5SJed Brown     The FULL factorization is
1357ab1df9f5SJed Brown 
1358ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1359ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1360ab1df9f5SJed Brown 
13616be4592eSBarry 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,
13626be4592eSBarry 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).
1363ab1df9f5SJed Brown 
13646be4592eSBarry 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
13656be4592eSBarry 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
13666be4592eSBarry 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,
13676be4592eSBarry 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.
1368ab1df9f5SJed Brown 
13696be4592eSBarry 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
13706be4592eSBarry 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).
1371ab1df9f5SJed Brown 
1372ab1df9f5SJed Brown     References:
1373ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1374ab1df9f5SJed Brown 
1375ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1376ab1df9f5SJed Brown 
1377ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1378ab1df9f5SJed Brown @*/
1379c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1380ab1df9f5SJed Brown {
1381ab1df9f5SJed Brown   PetscErrorCode ierr;
1382ab1df9f5SJed Brown 
1383ab1df9f5SJed Brown   PetscFunctionBegin;
1384ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1385c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1386ab1df9f5SJed Brown   PetscFunctionReturn(0);
1387ab1df9f5SJed Brown }
1388ab1df9f5SJed Brown 
1389ab1df9f5SJed Brown #undef __FUNCT__
1390c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1391c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1392ab1df9f5SJed Brown {
1393ab1df9f5SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1394ab1df9f5SJed Brown 
1395ab1df9f5SJed Brown   PetscFunctionBegin;
1396ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1397ab1df9f5SJed Brown   PetscFunctionReturn(0);
1398ab1df9f5SJed Brown }
1399ab1df9f5SJed Brown 
1400ab1df9f5SJed Brown #undef __FUNCT__
140130ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
140230ad9308SMatthew Knepley /*@C
14038c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
140430ad9308SMatthew Knepley 
140530ad9308SMatthew Knepley    Collective on KSP
140630ad9308SMatthew Knepley 
140730ad9308SMatthew Knepley    Input Parameter:
140830ad9308SMatthew Knepley .  pc - the preconditioner context
140930ad9308SMatthew Knepley 
141030ad9308SMatthew Knepley    Output Parameters:
1411a04f6461SBarry Smith +  A00 - the (0,0) block
1412a04f6461SBarry Smith .  A01 - the (0,1) block
1413a04f6461SBarry Smith .  A10 - the (1,0) block
1414a04f6461SBarry Smith -  A11 - the (1,1) block
141530ad9308SMatthew Knepley 
141630ad9308SMatthew Knepley    Level: advanced
141730ad9308SMatthew Knepley 
141830ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
141930ad9308SMatthew Knepley @*/
1420a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
142130ad9308SMatthew Knepley {
142230ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
142330ad9308SMatthew Knepley 
142430ad9308SMatthew Knepley   PetscFunctionBegin;
14250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1426c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1427a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1428a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1429a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1430a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
143130ad9308SMatthew Knepley   PetscFunctionReturn(0);
143230ad9308SMatthew Knepley }
143330ad9308SMatthew Knepley 
143479416396SBarry Smith EXTERN_C_BEGIN
143579416396SBarry Smith #undef __FUNCT__
143679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
14377087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
143879416396SBarry Smith {
143979416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1440e69d4d44SBarry Smith   PetscErrorCode ierr;
144179416396SBarry Smith 
144279416396SBarry Smith   PetscFunctionBegin;
144379416396SBarry Smith   jac->type = type;
14443b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
14453b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
14463b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1447e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1448e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1449c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1450e69d4d44SBarry Smith 
14513b224e63SBarry Smith   } else {
14523b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
14533b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1454e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14559e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1456c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
14573b224e63SBarry Smith   }
145879416396SBarry Smith   PetscFunctionReturn(0);
145979416396SBarry Smith }
146079416396SBarry Smith EXTERN_C_END
146179416396SBarry Smith 
146251f519a2SBarry Smith EXTERN_C_BEGIN
146351f519a2SBarry Smith #undef __FUNCT__
146451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
14657087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
146651f519a2SBarry Smith {
146751f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
146851f519a2SBarry Smith 
146951f519a2SBarry Smith   PetscFunctionBegin;
147065e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
147165e19b50SBarry 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);
147251f519a2SBarry Smith   jac->bs = bs;
147351f519a2SBarry Smith   PetscFunctionReturn(0);
147451f519a2SBarry Smith }
147551f519a2SBarry Smith EXTERN_C_END
147651f519a2SBarry Smith 
147779416396SBarry Smith #undef __FUNCT__
147879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1479bc08b0f1SBarry Smith /*@
148079416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
148179416396SBarry Smith 
148279416396SBarry Smith    Collective on PC
148379416396SBarry Smith 
148479416396SBarry Smith    Input Parameter:
148579416396SBarry Smith .  pc - the preconditioner context
148681540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
148779416396SBarry Smith 
148879416396SBarry Smith    Options Database Key:
1489a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
149079416396SBarry Smith 
149179416396SBarry Smith    Level: Developer
149279416396SBarry Smith 
149379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
149479416396SBarry Smith 
149579416396SBarry Smith .seealso: PCCompositeSetType()
149679416396SBarry Smith 
149779416396SBarry Smith @*/
14987087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
149979416396SBarry Smith {
15004ac538c5SBarry Smith   PetscErrorCode ierr;
150179416396SBarry Smith 
150279416396SBarry Smith   PetscFunctionBegin;
15030700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15044ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
150579416396SBarry Smith   PetscFunctionReturn(0);
150679416396SBarry Smith }
150779416396SBarry Smith 
15080971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
15090971522cSBarry Smith /*MC
1510a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1511a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
15120971522cSBarry Smith 
1513edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1514edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
15150971522cSBarry Smith 
1516a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
151769a612a9SBarry Smith          and set the options directly on the resulting KSP object
15180971522cSBarry Smith 
15190971522cSBarry Smith    Level: intermediate
15200971522cSBarry Smith 
152179416396SBarry Smith    Options Database Keys:
152281540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
152381540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
152481540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
152581540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
152681540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1527e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1528435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
152979416396SBarry Smith 
15305d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
15315d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
15325d4c12cdSJungho Lee 
1533c8a0d604SMatthew G Knepley    Notes:
1534c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1535d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1536d32f9abdSBarry Smith 
1537d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1538d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1539d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1540d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1541d32f9abdSBarry Smith 
1542c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1543c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1544c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1545c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1546c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1547a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1548c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1549c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1550c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1551a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1552c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1553c9c6ffaaSJed Brown      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
15545668aaf4SBarry Smith      diag gives
1555c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1556c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
15575668aaf4SBarry 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
1558c8a0d604SMatthew G Knepley $              (  A00   0 )
1559c8a0d604SMatthew G Knepley $              (  A10   S )
1560c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1561c8a0d604SMatthew G Knepley $              ( A00 A01 )
1562c8a0d604SMatthew G Knepley $              (  0   S  )
1563c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1564e69d4d44SBarry Smith 
1565edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1566edf189efSBarry Smith      is used automatically for a second block.
1567edf189efSBarry Smith 
1568ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1569ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1570ff218e97SBarry Smith 
1571ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1572ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1573ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
15740716a85fSBarry Smith 
1575a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
15760971522cSBarry Smith 
15777e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1578e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
15790971522cSBarry Smith M*/
15800971522cSBarry Smith 
15810971522cSBarry Smith EXTERN_C_BEGIN
15820971522cSBarry Smith #undef __FUNCT__
15830971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
15847087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
15850971522cSBarry Smith {
15860971522cSBarry Smith   PetscErrorCode ierr;
15870971522cSBarry Smith   PC_FieldSplit  *jac;
15880971522cSBarry Smith 
15890971522cSBarry Smith   PetscFunctionBegin;
159038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
15910971522cSBarry Smith   jac->bs        = -1;
15920971522cSBarry Smith   jac->nsplits   = 0;
15933e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1594e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1595c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
159651f519a2SBarry Smith 
15970971522cSBarry Smith   pc->data     = (void*)jac;
15980971522cSBarry Smith 
15990971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1600421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
16010971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1602574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
16030971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
16040971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
16050971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
16060971522cSBarry Smith   pc->ops->applyrichardson   = 0;
16070971522cSBarry Smith 
160869a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
160969a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
16100971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
16110971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1612704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1613704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
161479416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
161579416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
161651f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
161751f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
16180971522cSBarry Smith   PetscFunctionReturn(0);
16190971522cSBarry Smith }
16200971522cSBarry Smith EXTERN_C_END
16210971522cSBarry Smith 
16220971522cSBarry Smith 
1623a541d17aSBarry Smith 
1624