xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 20252d06f01a4399edd2208911327545c80253b1)
1dba47a55SKris Buschelman 
2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.h" I*/
40971522cSBarry Smith 
5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
7c5d2311dSJed Brown 
80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
90971522cSBarry Smith struct _PC_FieldSplitLink {
1069a612a9SBarry Smith   KSP               ksp;
110971522cSBarry Smith   Vec               x,y;
12db4c96c1SJed Brown   char              *splitname;
130971522cSBarry Smith   PetscInt          nfields;
145d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
151b9fc7fcSBarry Smith   VecScatter        sctx;
165d4c12cdSJungho Lee   IS                is,is_col;
1751f519a2SBarry Smith   PC_FieldSplitLink next,previous;
180971522cSBarry Smith };
190971522cSBarry Smith 
200971522cSBarry Smith typedef struct {
2168dd23aaSBarry Smith   PCCompositeType                    type;
22ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
23ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
24ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
2530ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
2630ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
2779416396SBarry Smith   Vec                                *x,*y,w1,w2;
28519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
29519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3030ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
31ace3abfcSBarry Smith   PetscBool                          issetup;
3230ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3330ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
3430ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
35a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
36084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
37084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
38c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
3930ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4097bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4163ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
42c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
430971522cSBarry Smith } PC_FieldSplit;
440971522cSBarry Smith 
4516913363SBarry Smith /*
4616913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
4716913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
4816913363SBarry Smith    PC you could change this.
4916913363SBarry Smith */
50084e4875SJed Brown 
51e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
52084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
53084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
54084e4875SJed Brown {
55084e4875SJed Brown   switch (jac->schurpre) {
56084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
57084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
58084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
59084e4875SJed Brown     default:
60084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
61084e4875SJed Brown   }
62084e4875SJed Brown }
63084e4875SJed Brown 
64084e4875SJed Brown 
650971522cSBarry Smith #undef __FUNCT__
660971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
670971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
680971522cSBarry Smith {
690971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
700971522cSBarry Smith   PetscErrorCode    ierr;
71ace3abfcSBarry Smith   PetscBool         iascii;
720971522cSBarry Smith   PetscInt          i,j;
735a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
740971522cSBarry Smith 
750971522cSBarry Smith   PetscFunctionBegin;
76251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
770971522cSBarry Smith   if (iascii) {
782c9966d7SBarry Smith     if (jac->bs > 0) {
7951f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
802c9966d7SBarry Smith     } else {
812c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
822c9966d7SBarry Smith     }
83a3df900dSMatthew G Knepley     if (jac->realdiagonal) {
84a3df900dSMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr);
85a3df900dSMatthew G Knepley     }
8669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
870971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
880971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
891ab39975SBarry Smith       if (ilink->fields) {
900971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9179416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
925a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9379416396SBarry Smith 	  if (j > 0) {
9479416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9579416396SBarry Smith 	  }
965a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
970971522cSBarry Smith 	}
980971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1001ab39975SBarry Smith       } else {
1011ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1021ab39975SBarry Smith       }
1035a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1045a9f2f41SSatish Balay       ilink = ilink->next;
1050971522cSBarry Smith     }
1060971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1070971522cSBarry Smith   } else {
10865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1090971522cSBarry Smith   }
1100971522cSBarry Smith   PetscFunctionReturn(0);
1110971522cSBarry Smith }
1120971522cSBarry Smith 
1130971522cSBarry Smith #undef __FUNCT__
1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1163b224e63SBarry Smith {
1173b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1183b224e63SBarry Smith   PetscErrorCode    ierr;
119ace3abfcSBarry Smith   PetscBool         iascii;
1203b224e63SBarry Smith   PetscInt          i,j;
1213b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1223b224e63SBarry Smith   KSP               ksp;
1233b224e63SBarry Smith 
1243b224e63SBarry Smith   PetscFunctionBegin;
125251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1263b224e63SBarry Smith   if (iascii) {
1272c9966d7SBarry Smith     if (jac->bs > 0) {
128c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1292c9966d7SBarry Smith     } else {
1302c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1312c9966d7SBarry Smith     }
132a3df900dSMatthew G Knepley     if (jac->realdiagonal) {
133a3df900dSMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr);
134a3df900dSMatthew G Knepley     }
1353e8b8b31SMatthew G Knepley     switch(jac->schurpre) {
1363e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
1373e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break;
138b02e2d75SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_DIAG:
1393e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break;
1403e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1413e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1423e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1433e8b8b31SMatthew G Knepley       } else {
1443e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);
1453e8b8b31SMatthew G Knepley       }
1463e8b8b31SMatthew G Knepley       break;
1473e8b8b31SMatthew G Knepley     default:
1483e8b8b31SMatthew G Knepley       SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1493e8b8b31SMatthew G Knepley     }
1503b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1513b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1523b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1533b224e63SBarry Smith       if (ilink->fields) {
1543b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1553b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1563b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1573b224e63SBarry Smith 	  if (j > 0) {
1583b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1593b224e63SBarry Smith 	  }
1603b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1613b224e63SBarry Smith 	}
1623b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1633b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1643b224e63SBarry Smith       } else {
1653b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1663b224e63SBarry Smith       }
1673b224e63SBarry Smith       ilink = ilink->next;
1683b224e63SBarry Smith     }
169435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1703b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
17112cae6f2SJed Brown     if (jac->schur) {
17212cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1733b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
17412cae6f2SJed Brown     } else {
17512cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
17612cae6f2SJed Brown     }
1773b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
178435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1793b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18012cae6f2SJed Brown     if (jac->kspschur) {
1813b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
18212cae6f2SJed Brown     } else {
18312cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
18412cae6f2SJed Brown     }
1853b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1863b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1873b224e63SBarry Smith   } else {
18865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1893b224e63SBarry Smith   }
1903b224e63SBarry Smith   PetscFunctionReturn(0);
1913b224e63SBarry Smith }
1923b224e63SBarry Smith 
1933b224e63SBarry Smith #undef __FUNCT__
1946c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1956c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1966c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1976c924f48SJed Brown {
1986c924f48SJed Brown   PetscErrorCode ierr;
1996c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2005d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
2015d4c12cdSJungho Lee   PetscBool      flg,flg_col;
2025d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
2036c924f48SJed Brown 
2046c924f48SJed Brown   PetscFunctionBegin;
2056c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
2065d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
2076c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
2086c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2096c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
2105d4c12cdSJungho Lee     ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
2116c924f48SJed Brown     nfields = jac->bs;
21229499fbbSJungho Lee     nfields_col = jac->bs;
2136c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
2145d4c12cdSJungho Lee     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
2156c924f48SJed Brown     if (!flg) break;
2165d4c12cdSJungho Lee     else if (flg && !flg_col) {
2175d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2185d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
2195d4c12cdSJungho Lee     }
2205d4c12cdSJungho Lee     else {
2215d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2225d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
2235d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
2245d4c12cdSJungho Lee     }
2256c924f48SJed Brown   }
2266c924f48SJed Brown   if (i > 0) {
2276c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
2286c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2296c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2306c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2316c924f48SJed Brown   }
2326c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2335d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2346c924f48SJed Brown   PetscFunctionReturn(0);
2356c924f48SJed Brown }
2366c924f48SJed Brown 
2376c924f48SJed Brown #undef __FUNCT__
23869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
23969a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2400971522cSBarry Smith {
2410971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2420971522cSBarry Smith   PetscErrorCode    ierr;
2435a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2447287d2a3SDmitry Karpeev   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE;
2456c924f48SJed Brown   PetscInt          i;
2460971522cSBarry Smith 
2470971522cSBarry Smith   PetscFunctionBegin;
2487287d2a3SDmitry Karpeev   /*
2497287d2a3SDmitry Karpeev    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
2507287d2a3SDmitry Karpeev    Should probably be rewritten.
2517287d2a3SDmitry Karpeev    */
2527287d2a3SDmitry Karpeev   if (!ilink || jac->reset) {
253d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
254d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
255bafc1b83SMatthew G Knepley       PetscInt     numFields, f, i, j;
2560784a22cSJed Brown       char       **fieldNames;
2577b62db95SJungho Lee       IS          *fields;
258e7c4fc90SDmitry Karpeev       DM          *dms;
259bafc1b83SMatthew G Knepley       DM           subdm[128];
260bafc1b83SMatthew G Knepley       PetscBool    flg;
261bafc1b83SMatthew G Knepley 
26216621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
263bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
264bafc1b83SMatthew G Knepley       for(i = 0, flg = PETSC_TRUE; ; i++) {
265bafc1b83SMatthew G Knepley         PetscInt ifields[128];
266bafc1b83SMatthew G Knepley         IS       compField;
267bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
268bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
269bafc1b83SMatthew G Knepley 
270bafc1b83SMatthew G Knepley         ierr = PetscSNPrintf(optionname, sizeof optionname, "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
271bafc1b83SMatthew G Knepley         ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
272bafc1b83SMatthew G Knepley         if (!flg) break;
273bafc1b83SMatthew G Knepley         if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
274bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
275bafc1b83SMatthew G Knepley         if (nfields == 1) {
276bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
2774d2389d7SMatthew G Knepley           /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr);
2784d2389d7SMatthew G Knepley              ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */
279bafc1b83SMatthew G Knepley         } else {
280bafc1b83SMatthew G Knepley           ierr = PetscSNPrintf(splitname, sizeof splitname, "%D", i);CHKERRQ(ierr);
281bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
2824d2389d7SMatthew G Knepley           /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr);
2834d2389d7SMatthew G Knepley              ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */
2847287d2a3SDmitry Karpeev         }
285bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
286bafc1b83SMatthew G Knepley         for(j = 0; j < nfields; ++j) {
287bafc1b83SMatthew G Knepley           f = ifields[j];
2887b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
2897b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
2907b62db95SJungho Lee         }
291bafc1b83SMatthew G Knepley       }
292bafc1b83SMatthew G Knepley       if (i == 0) {
293bafc1b83SMatthew G Knepley         for(f = 0; f < numFields; ++f) {
294bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
295bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
296bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
297bafc1b83SMatthew G Knepley         }
298bafc1b83SMatthew G Knepley       } else {
299bafc1b83SMatthew G Knepley         ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr);
300bafc1b83SMatthew G Knepley         for(j = 0; j < i; ++j) {
301bafc1b83SMatthew G Knepley           dms[j] = subdm[j];
302bafc1b83SMatthew G Knepley         }
303bafc1b83SMatthew G Knepley       }
3047b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
3057b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
306e7c4fc90SDmitry Karpeev       if (dms) {
3078b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
308bafc1b83SMatthew G Knepley         for(ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
3097287d2a3SDmitry Karpeev           const char *prefix;
3107287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr);
3117287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);     CHKERRQ(ierr);
3127b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
3137b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
3147287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr);
315e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]); CHKERRQ(ierr);
3162fa5ba8aSJed Brown         }
3177b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
3188b8307b2SJed Brown       }
31966ffff09SJed Brown     } else {
320521d7252SBarry Smith       if (jac->bs <= 0) {
321704ba839SBarry Smith         if (pc->pmat) {
322521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
323704ba839SBarry Smith         } else {
324704ba839SBarry Smith           jac->bs = 1;
325704ba839SBarry Smith         }
326521d7252SBarry Smith       }
327d32f9abdSBarry Smith 
3287287d2a3SDmitry Karpeev       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr);
3296ce1633cSBarry Smith       if (stokes) {
3306ce1633cSBarry Smith         IS       zerodiags,rest;
3316ce1633cSBarry Smith         PetscInt nmin,nmax;
3326ce1633cSBarry Smith 
3336ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
3346ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
3356ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
3367287d2a3SDmitry Karpeev         if(jac->reset) {
3377287d2a3SDmitry Karpeev           jac->head->is       = rest;
3387287d2a3SDmitry Karpeev           jac->head->next->is = zerodiags;
3397287d2a3SDmitry Karpeev         }
3407287d2a3SDmitry Karpeev         else {
3416ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
3426ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
3437287d2a3SDmitry Karpeev         }
344fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
345fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
3466ce1633cSBarry Smith       } else {
3477287d2a3SDmitry Karpeev         if(jac->reset)
3487287d2a3SDmitry Karpeev           SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
3497287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
350d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
351d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
3526c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
3536c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
354d32f9abdSBarry Smith         }
3557287d2a3SDmitry Karpeev         if (fieldsplit_default || !jac->splitdefined) {
356d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
357db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
3586c924f48SJed Brown             char splitname[8];
3596c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
3605d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
36179416396SBarry Smith           }
3625d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
363521d7252SBarry Smith         }
36466ffff09SJed Brown       }
3656ce1633cSBarry Smith     }
366edf189efSBarry Smith   } else if (jac->nsplits == 1) {
367edf189efSBarry Smith     if (ilink->is) {
368edf189efSBarry Smith       IS       is2;
369edf189efSBarry Smith       PetscInt nmin,nmax;
370edf189efSBarry Smith 
371edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
372edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
373db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
374fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
3757b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
376edf189efSBarry Smith   }
377d0af7cd3SBarry Smith 
378d0af7cd3SBarry Smith 
3797b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
38069a612a9SBarry Smith   PetscFunctionReturn(0);
38169a612a9SBarry Smith }
38269a612a9SBarry Smith 
383514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg);
384514bf10dSMatthew G Knepley 
38569a612a9SBarry Smith #undef __FUNCT__
38669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
38769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
38869a612a9SBarry Smith {
38969a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
39069a612a9SBarry Smith   PetscErrorCode    ierr;
3915a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
3922c9966d7SBarry Smith   PetscInt          i,nsplit;
39369a612a9SBarry Smith   MatStructure      flag = pc->flag;
3945f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
39569a612a9SBarry Smith 
39669a612a9SBarry Smith   PetscFunctionBegin;
39769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
39897bbdb24SBarry Smith   nsplit = jac->nsplits;
3995a9f2f41SSatish Balay   ilink  = jac->head;
40097bbdb24SBarry Smith 
40197bbdb24SBarry Smith   /* get the matrices for each split */
402704ba839SBarry Smith   if (!jac->issetup) {
4031b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
40497bbdb24SBarry Smith 
405704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
406704ba839SBarry Smith 
4075d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4082c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4092c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
4102c9966d7SBarry Smith     }
41151f519a2SBarry Smith     bs     = jac->bs;
41297bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
4131b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4141b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4151b9fc7fcSBarry Smith       if (jac->defaultsplit) {
416704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
4175f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
418704ba839SBarry Smith       } else if (!ilink->is) {
419ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4205f4ab4e1SJungho Lee           PetscInt   *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
4215a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
4225f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
4231b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
4241b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
4251b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
4265f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
42797bbdb24SBarry Smith             }
42897bbdb24SBarry Smith           }
429d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
4305f4ab4e1SJungho Lee           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
431ccb205f8SBarry Smith         } else {
432704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
4335f4ab4e1SJungho Lee           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
434ccb205f8SBarry Smith        }
4353e197d65SBarry Smith       }
436704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
4375f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
4385f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
439704ba839SBarry Smith       ilink = ilink->next;
4401b9fc7fcSBarry Smith     }
4411b9fc7fcSBarry Smith   }
4421b9fc7fcSBarry Smith 
443704ba839SBarry Smith   ilink  = jac->head;
44497bbdb24SBarry Smith   if (!jac->pmat) {
445bdddcaaaSMatthew G Knepley     Vec xtmp;
446bdddcaaaSMatthew G Knepley 
447bdddcaaaSMatthew G Knepley     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
448cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
449bdddcaaaSMatthew G Knepley     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
450cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
451bdddcaaaSMatthew G Knepley       MatNullSpace sp;
452bdddcaaaSMatthew G Knepley 
453a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
454a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr);
455a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
456a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
457a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
458a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
459a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
460a3df900dSMatthew G Knepley         }
461a3df900dSMatthew G Knepley       } else {
4625f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
463a3df900dSMatthew G Knepley       }
464bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
465bdddcaaaSMatthew G Knepley       ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
466bdddcaaaSMatthew G Knepley       ilink->x = jac->x[i]; ilink->y = jac->y[i];
467bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
468bdddcaaaSMatthew G Knepley       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
469bdddcaaaSMatthew G Knepley       /* HACK: Check for the constant null space */
470bdddcaaaSMatthew G Knepley       ierr = MatGetNullSpace(pc->pmat, &sp);CHKERRQ(ierr);
471bdddcaaaSMatthew G Knepley       if (sp) {
472bdddcaaaSMatthew G Knepley         MatNullSpace subsp;
473bdddcaaaSMatthew G Knepley         Vec          ftmp, gtmp;
4749583d628SBarry Smith         PetscReal    norm;
475bdddcaaaSMatthew G Knepley         PetscInt     N;
476bdddcaaaSMatthew G Knepley 
477bdddcaaaSMatthew G Knepley         ierr = MatGetVecs(pc->pmat,     &gtmp, PETSC_NULL);CHKERRQ(ierr);
478bdddcaaaSMatthew G Knepley         ierr = MatGetVecs(jac->pmat[i], &ftmp, PETSC_NULL);CHKERRQ(ierr);
479bdddcaaaSMatthew G Knepley         ierr = VecGetSize(ftmp, &N);CHKERRQ(ierr);
480bdddcaaaSMatthew G Knepley         ierr = VecSet(ftmp, 1.0/N);CHKERRQ(ierr);
481bdddcaaaSMatthew G Knepley         ierr = VecScatterBegin(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
482bdddcaaaSMatthew G Knepley         ierr = VecScatterEnd(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
483bdddcaaaSMatthew G Knepley         ierr = MatNullSpaceRemove(sp, gtmp, PETSC_NULL);CHKERRQ(ierr);
484bdddcaaaSMatthew G Knepley         ierr = VecNorm(gtmp, NORM_2, &norm);CHKERRQ(ierr);
485bdddcaaaSMatthew G Knepley         if (norm < 1.0e-10) {
486bdddcaaaSMatthew G Knepley           ierr  = MatNullSpaceCreate(((PetscObject)pc)->comm, PETSC_TRUE, 0, PETSC_NULL, &subsp);CHKERRQ(ierr);
487bdddcaaaSMatthew G Knepley           ierr  = MatSetNullSpace(jac->pmat[i], subsp);CHKERRQ(ierr);
488bdddcaaaSMatthew G Knepley           ierr  = MatNullSpaceDestroy(&subsp);CHKERRQ(ierr);
489bdddcaaaSMatthew G Knepley         }
490bdddcaaaSMatthew G Knepley         ierr = VecDestroy(&ftmp);CHKERRQ(ierr);
491bdddcaaaSMatthew G Knepley         ierr = VecDestroy(&gtmp);CHKERRQ(ierr);
492bdddcaaaSMatthew G Knepley       }
493ed1f0337SMatthew G Knepley       /* Check for null space attached to IS */
494bafc1b83SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr);
495bafc1b83SMatthew G Knepley       if (sp) {
496bafc1b83SMatthew G Knepley         ierr  = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
497bafc1b83SMatthew G Knepley       }
498ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr);
499ed1f0337SMatthew G Knepley       if (sp) {
500ed1f0337SMatthew G Knepley         ierr  = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
501ed1f0337SMatthew G Knepley       }
502704ba839SBarry Smith       ilink = ilink->next;
503cf502942SBarry Smith     }
504bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
50597bbdb24SBarry Smith   } else {
506cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
507a3df900dSMatthew G Knepley       Mat pmat;
508a3df900dSMatthew G Knepley 
509a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
510a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr);
511a3df900dSMatthew G Knepley       if (!pmat) {
5125f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
513a3df900dSMatthew G Knepley       }
514704ba839SBarry Smith       ilink = ilink->next;
515cf502942SBarry Smith     }
51697bbdb24SBarry Smith   }
517519d70e2SJed Brown   if (jac->realdiagonal) {
518519d70e2SJed Brown     ilink = jac->head;
519519d70e2SJed Brown     if (!jac->mat) {
520519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
521519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5225f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
523519d70e2SJed Brown         ilink = ilink->next;
524519d70e2SJed Brown       }
525519d70e2SJed Brown     } else {
526519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5275f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
528519d70e2SJed Brown         ilink = ilink->next;
529519d70e2SJed Brown       }
530519d70e2SJed Brown     }
531519d70e2SJed Brown   } else {
532519d70e2SJed Brown     jac->mat = jac->pmat;
533519d70e2SJed Brown   }
53497bbdb24SBarry Smith 
5356c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
53668dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
53768dd23aaSBarry Smith     ilink  = jac->head;
53868dd23aaSBarry Smith     if (!jac->Afield) {
53968dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
54068dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5414aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
54268dd23aaSBarry Smith         ilink = ilink->next;
54368dd23aaSBarry Smith       }
54468dd23aaSBarry Smith     } else {
54568dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5464aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
54768dd23aaSBarry Smith         ilink = ilink->next;
54868dd23aaSBarry Smith       }
54968dd23aaSBarry Smith     }
55068dd23aaSBarry Smith   }
55168dd23aaSBarry Smith 
5523b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
5533b224e63SBarry Smith     IS       ccis;
5544aa3045dSJed Brown     PetscInt rstart,rend;
555093c86ffSJed Brown     char     lscname[256];
556093c86ffSJed Brown     PetscObject LSC_L;
557e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
55868dd23aaSBarry Smith 
559e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
560e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
561e6cab6aaSJed Brown 
5623b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
5633b224e63SBarry Smith     if (jac->schur) {
5643b224e63SBarry Smith       ilink = jac->head;
56549bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5664aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
567fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5683b224e63SBarry Smith       ilink = ilink->next;
56949bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5704aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
571fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
572a3df900dSMatthew G Knepley       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr);
573084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
5743b224e63SBarry Smith      } else {
5751cee3971SBarry Smith       KSP ksp;
576bafc1b83SMatthew G Knepley       const char  *Dprefix;
577bafc1b83SMatthew G Knepley       char schurprefix[256];
578514bf10dSMatthew G Knepley       char schurtestoption[256];
579bdddcaaaSMatthew G Knepley       MatNullSpace sp;
580514bf10dSMatthew G Knepley       PetscBool    flg;
5813b224e63SBarry Smith 
582a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
5833b224e63SBarry Smith       ilink = jac->head;
58449bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5854aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
586fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5873b224e63SBarry Smith       ilink = ilink->next;
58849bb4cd7SJungho Lee       ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
5894aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
590fcfd50ebSBarry Smith       ierr = ISDestroy(&ccis);CHKERRQ(ierr);
591*20252d06SBarry Smith 
592*20252d06SBarry Smith       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */
593*20252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
594*20252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
595*20252d06SBarry Smith       ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr);
596*20252d06SBarry Smith       ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
597*20252d06SBarry Smith       /* Indent this deeper to emphasize the "inner" nature of this solver. */
598*20252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr);
599*20252d06SBarry Smith       ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr);
600*20252d06SBarry Smith       ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
601*20252d06SBarry Smith 
602bdddcaaaSMatthew G Knepley       ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr);
603*20252d06SBarry Smith       if (sp) {
604*20252d06SBarry Smith         ierr  = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
605*20252d06SBarry Smith       }
606*20252d06SBarry Smith 
607*20252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
608514bf10dSMatthew G Knepley       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr);
609514bf10dSMatthew G Knepley       if (flg) {
610514bf10dSMatthew G Knepley         DM dmInner;
611514bf10dSMatthew G Knepley 
612514bf10dSMatthew G Knepley         /* Set DM for new solver */
613bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
614bafc1b83SMatthew G Knepley         ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr);
6157287d2a3SDmitry Karpeev         ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
616514bf10dSMatthew G Knepley       } else {
617514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr);
618514bf10dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
619bafc1b83SMatthew G Knepley       }
62020b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
6213b224e63SBarry Smith 
6223b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);               CHKERRQ(ierr);
6239005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
624*20252d06SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);           CHKERRQ(ierr);
625084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
626084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
627084e4875SJed Brown         PC pc;
628cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
629084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
630084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
631e69d4d44SBarry Smith       }
6327287d2a3SDmitry Karpeev       ierr  = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr);
6337287d2a3SDmitry Karpeev       ierr  = KSPSetOptionsPrefix(jac->kspschur,         Dprefix); CHKERRQ(ierr);
6343b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
63520b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
63620b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
6373b224e63SBarry Smith     }
638093c86ffSJed Brown 
639093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
640093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
641093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
642093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
643093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
644093c86ffSJed Brown     ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
645093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
646093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
647093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
6483b224e63SBarry Smith   } else {
64997bbdb24SBarry Smith     /* set up the individual PCs */
65097bbdb24SBarry Smith     i    = 0;
6515a9f2f41SSatish Balay     ilink = jac->head;
6525a9f2f41SSatish Balay     while (ilink) {
653519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
6543b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
655c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
6565a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
65797bbdb24SBarry Smith       i++;
6585a9f2f41SSatish Balay       ilink = ilink->next;
6590971522cSBarry Smith     }
6603b224e63SBarry Smith   }
6613b224e63SBarry Smith 
662c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
6630971522cSBarry Smith   PetscFunctionReturn(0);
6640971522cSBarry Smith }
6650971522cSBarry Smith 
6665a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
667ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
668ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
6695a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
670ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
671ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
67279416396SBarry Smith 
6730971522cSBarry Smith #undef __FUNCT__
6743b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
6753b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
6763b224e63SBarry Smith {
6773b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6783b224e63SBarry Smith   PetscErrorCode    ierr;
6793b224e63SBarry Smith   KSP               ksp;
6803b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
6813b224e63SBarry Smith 
6823b224e63SBarry Smith   PetscFunctionBegin;
6833b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
6843b224e63SBarry Smith 
685c5d2311dSJed Brown   switch (jac->schurfactorization) {
686c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
687a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
688c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
689c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
690c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
691c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
692c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
693c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
694c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
695c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
696c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
697c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
698c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
699c5d2311dSJed Brown       break;
700c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
701a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
702c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
703c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
704c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
705c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
706c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
707c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
708c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
709c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
710c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
711c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
712c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
713c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
714c5d2311dSJed Brown       break;
715c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
716a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
717c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
718c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
719c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
720c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
721c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
722c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
723c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
724c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
725c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
726c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
727c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
728c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
729c5d2311dSJed Brown       break;
730c9c6ffaaSJed Brown     case PC_FIELDSPLIT_SCHUR_FACT_FULL:
731a04f6461SBarry 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 */
7323b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7333b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7343b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
7353b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
7363b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
7373b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7383b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7393b224e63SBarry Smith 
7403b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
7413b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7423b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7433b224e63SBarry Smith 
7443b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
7453b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
7463b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
7473b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7483b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
749c5d2311dSJed Brown   }
7503b224e63SBarry Smith   PetscFunctionReturn(0);
7513b224e63SBarry Smith }
7523b224e63SBarry Smith 
7533b224e63SBarry Smith #undef __FUNCT__
7540971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
7550971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
7560971522cSBarry Smith {
7570971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7580971522cSBarry Smith   PetscErrorCode    ierr;
7595a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
760939b8a20SBarry Smith   PetscInt          cnt,bs;
7610971522cSBarry Smith 
7620971522cSBarry Smith   PetscFunctionBegin;
7634442daceSBarry Smith   x->map->bs = jac->bs;
7644442daceSBarry Smith   y->map->bs = jac->bs;
76551f519a2SBarry Smith   CHKMEMQ;
76679416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
7671b9fc7fcSBarry Smith     if (jac->defaultsplit) {
768939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
769939b8a20SBarry 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);
770939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
771939b8a20SBarry 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);
7720971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
7735a9f2f41SSatish Balay       while (ilink) {
7745a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7755a9f2f41SSatish Balay         ilink = ilink->next;
7760971522cSBarry Smith       }
7770971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
7781b9fc7fcSBarry Smith     } else {
779efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
7805a9f2f41SSatish Balay       while (ilink) {
7815a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7825a9f2f41SSatish Balay         ilink = ilink->next;
7831b9fc7fcSBarry Smith       }
7841b9fc7fcSBarry Smith     }
78516913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
78679416396SBarry Smith     if (!jac->w1) {
78779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
78879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
78979416396SBarry Smith     }
790efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
7915a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7923e197d65SBarry Smith     cnt = 1;
7935a9f2f41SSatish Balay     while (ilink->next) {
7945a9f2f41SSatish Balay       ilink = ilink->next;
7953e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7963e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7973e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7983e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7993e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8003e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
8013e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8023e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8033e197d65SBarry Smith     }
80451f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
80511755939SBarry Smith       cnt -= 2;
80651f519a2SBarry Smith       while (ilink->previous) {
80751f519a2SBarry Smith         ilink = ilink->previous;
80811755939SBarry Smith         /* compute the residual only over the part of the vector needed */
80911755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
81011755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
81111755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81211755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81311755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
81411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
81511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
81651f519a2SBarry Smith       }
81711755939SBarry Smith     }
81865e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
81951f519a2SBarry Smith   CHKMEMQ;
8200971522cSBarry Smith   PetscFunctionReturn(0);
8210971522cSBarry Smith }
8220971522cSBarry Smith 
823421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
824ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
825ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
826421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
827ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
828ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
829421e10b8SBarry Smith 
830421e10b8SBarry Smith #undef __FUNCT__
8318c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
832421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
833421e10b8SBarry Smith {
834421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
835421e10b8SBarry Smith   PetscErrorCode    ierr;
836421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
837939b8a20SBarry Smith   PetscInt          bs;
838421e10b8SBarry Smith 
839421e10b8SBarry Smith   PetscFunctionBegin;
840421e10b8SBarry Smith   CHKMEMQ;
841421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
842421e10b8SBarry Smith     if (jac->defaultsplit) {
843939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
844939b8a20SBarry 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);
845939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
846939b8a20SBarry 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);
847421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
848421e10b8SBarry Smith       while (ilink) {
849421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
850421e10b8SBarry Smith 	ilink = ilink->next;
851421e10b8SBarry Smith       }
852421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
853421e10b8SBarry Smith     } else {
854421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
855421e10b8SBarry Smith       while (ilink) {
856421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
857421e10b8SBarry Smith 	ilink = ilink->next;
858421e10b8SBarry Smith       }
859421e10b8SBarry Smith     }
860421e10b8SBarry Smith   } else {
861421e10b8SBarry Smith     if (!jac->w1) {
862421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
863421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
864421e10b8SBarry Smith     }
865421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
866421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
867421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
868421e10b8SBarry Smith       while (ilink->next) {
869421e10b8SBarry Smith         ilink = ilink->next;
8709989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
871421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
872421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
873421e10b8SBarry Smith       }
874421e10b8SBarry Smith       while (ilink->previous) {
875421e10b8SBarry Smith         ilink = ilink->previous;
8769989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
877421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
878421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
879421e10b8SBarry Smith       }
880421e10b8SBarry Smith     } else {
881421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
882421e10b8SBarry Smith 	ilink = ilink->next;
883421e10b8SBarry Smith       }
884421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
885421e10b8SBarry Smith       while (ilink->previous) {
886421e10b8SBarry Smith 	ilink = ilink->previous;
8879989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
888421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
889421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
890421e10b8SBarry Smith       }
891421e10b8SBarry Smith     }
892421e10b8SBarry Smith   }
893421e10b8SBarry Smith   CHKMEMQ;
894421e10b8SBarry Smith   PetscFunctionReturn(0);
895421e10b8SBarry Smith }
896421e10b8SBarry Smith 
8970971522cSBarry Smith #undef __FUNCT__
898574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
899574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
9000971522cSBarry Smith {
9010971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9020971522cSBarry Smith   PetscErrorCode    ierr;
9035a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
9040971522cSBarry Smith 
9050971522cSBarry Smith   PetscFunctionBegin;
9065a9f2f41SSatish Balay   while (ilink) {
907574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
908fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
909fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
910fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
911fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
912b5787286SJed Brown     ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
9135a9f2f41SSatish Balay     next = ilink->next;
9145a9f2f41SSatish Balay     ilink = next;
9150971522cSBarry Smith   }
91605b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
917574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
918574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
919574deadeSBarry Smith   } else if (jac->mat) {
920574deadeSBarry Smith     jac->mat = PETSC_NULL;
921574deadeSBarry Smith   }
92297bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
92368dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
9246bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
9256bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
9266bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
9276bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
9286bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
9296bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
9306bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
93163ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
932574deadeSBarry Smith   PetscFunctionReturn(0);
933574deadeSBarry Smith }
934574deadeSBarry Smith 
935574deadeSBarry Smith #undef __FUNCT__
936574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
937574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
938574deadeSBarry Smith {
939574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
940574deadeSBarry Smith   PetscErrorCode    ierr;
941574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
942574deadeSBarry Smith 
943574deadeSBarry Smith   PetscFunctionBegin;
944574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
945574deadeSBarry Smith   while (ilink) {
9466bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
947574deadeSBarry Smith     next = ilink->next;
948574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
949574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
9505d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
951574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
952574deadeSBarry Smith     ilink = next;
953574deadeSBarry Smith   }
954574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
955c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
9567b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
9577b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
9587b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
9597b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
9607b62db95SJungho Lee   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
961ab1df9f5SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr);
962c9c6ffaaSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr);
9630971522cSBarry Smith   PetscFunctionReturn(0);
9640971522cSBarry Smith }
9650971522cSBarry Smith 
9660971522cSBarry Smith #undef __FUNCT__
9670971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
9680971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
9690971522cSBarry Smith {
9701b9fc7fcSBarry Smith   PetscErrorCode  ierr;
9716c924f48SJed Brown   PetscInt        bs;
972bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
9739dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
9743b224e63SBarry Smith   PCCompositeType ctype;
975e7c4fc90SDmitry Karpeev   DM              ddm;
976e7c4fc90SDmitry Karpeev   char            ddm_name[1025];
9771b9fc7fcSBarry Smith 
9780971522cSBarry Smith   PetscFunctionBegin;
9791b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
980e7c4fc90SDmitry Karpeev   if(pc->dm) {
981e7c4fc90SDmitry Karpeev     /* Allow the user to request a decomposition DM by name */
982e7c4fc90SDmitry Karpeev     ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
983731c8d9eSDmitry Karpeev     ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
984e7c4fc90SDmitry Karpeev     if(flg) {
98516621825SDmitry Karpeev       ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
986e7c4fc90SDmitry Karpeev       if(!ddm) {
98716621825SDmitry Karpeev         SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name);
988e7c4fc90SDmitry Karpeev       }
98916621825SDmitry Karpeev       ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr);
990e7c4fc90SDmitry Karpeev       ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
991e7c4fc90SDmitry Karpeev     }
992e7c4fc90SDmitry Karpeev   }
993acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
99451f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
99551f519a2SBarry Smith   if (flg) {
99651f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
99751f519a2SBarry Smith   }
998704ba839SBarry Smith 
999435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
1000c0adfefeSBarry Smith   if (stokes) {
1001c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1002c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1003c0adfefeSBarry Smith   }
1004c0adfefeSBarry Smith 
10053b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
10063b224e63SBarry Smith   if (flg) {
10073b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
10083b224e63SBarry Smith   }
1009c30613efSMatthew Knepley   /* Only setup fields once */
1010c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1011d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1012d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
10136c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
10146c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1015d32f9abdSBarry Smith   }
1016c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1017c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1018c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
1019c9c6ffaaSJed 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);
1020084e4875SJed 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);
1021c5d2311dSJed Brown   }
10221b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10230971522cSBarry Smith   PetscFunctionReturn(0);
10240971522cSBarry Smith }
10250971522cSBarry Smith 
10260971522cSBarry Smith /*------------------------------------------------------------------------------------*/
10270971522cSBarry Smith 
10280971522cSBarry Smith EXTERN_C_BEGIN
10290971522cSBarry Smith #undef __FUNCT__
10300971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
10315d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
10320971522cSBarry Smith {
103397bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10340971522cSBarry Smith   PetscErrorCode    ierr;
10355a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
103669a612a9SBarry Smith   char              prefix[128];
10375d4c12cdSJungho Lee   PetscInt          i;
10380971522cSBarry Smith 
10390971522cSBarry 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   }
104451f519a2SBarry Smith   for (i=0; i<n; i++) {
1045e32f2f54SBarry 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);
1046e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
104751f519a2SBarry Smith   }
1048704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1049a04f6461SBarry Smith   if (splitname) {
1050db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1051a04f6461SBarry Smith   } else {
1052a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1053a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1054a04f6461SBarry Smith   }
1055704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
10565a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
10575d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
10585d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
10595a9f2f41SSatish Balay   ilink->nfields = n;
10605a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
10617adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
1062*20252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
10635a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10649005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
106569a612a9SBarry Smith 
1066a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
10675a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
10680971522cSBarry Smith 
10690971522cSBarry Smith   if (!next) {
10705a9f2f41SSatish Balay     jac->head       = ilink;
107151f519a2SBarry Smith     ilink->previous = PETSC_NULL;
10720971522cSBarry Smith   } else {
10730971522cSBarry Smith     while (next->next) {
10740971522cSBarry Smith       next = next->next;
10750971522cSBarry Smith     }
10765a9f2f41SSatish Balay     next->next      = ilink;
107751f519a2SBarry Smith     ilink->previous = next;
10780971522cSBarry Smith   }
10790971522cSBarry Smith   jac->nsplits++;
10800971522cSBarry Smith   PetscFunctionReturn(0);
10810971522cSBarry Smith }
10820971522cSBarry Smith EXTERN_C_END
10830971522cSBarry Smith 
1084e69d4d44SBarry Smith EXTERN_C_BEGIN
1085e69d4d44SBarry Smith #undef __FUNCT__
1086e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
10877087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1088e69d4d44SBarry Smith {
1089e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1090e69d4d44SBarry Smith   PetscErrorCode ierr;
1091e69d4d44SBarry Smith 
1092e69d4d44SBarry Smith   PetscFunctionBegin;
1093519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1094e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1095e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
109613e0d083SBarry Smith   if (n) *n = jac->nsplits;
1097e69d4d44SBarry Smith   PetscFunctionReturn(0);
1098e69d4d44SBarry Smith }
1099e69d4d44SBarry Smith EXTERN_C_END
11000971522cSBarry Smith 
11010971522cSBarry Smith EXTERN_C_BEGIN
11020971522cSBarry Smith #undef __FUNCT__
110369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
11047087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
11050971522cSBarry Smith {
11060971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11070971522cSBarry Smith   PetscErrorCode    ierr;
11080971522cSBarry Smith   PetscInt          cnt = 0;
11095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
11100971522cSBarry Smith 
11110971522cSBarry Smith   PetscFunctionBegin;
11125d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
11135a9f2f41SSatish Balay   while (ilink) {
11145a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
11155a9f2f41SSatish Balay     ilink = ilink->next;
11160971522cSBarry Smith   }
11175d480477SMatthew 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);
111813e0d083SBarry Smith   if (n) *n = jac->nsplits;
11190971522cSBarry Smith   PetscFunctionReturn(0);
11200971522cSBarry Smith }
11210971522cSBarry Smith EXTERN_C_END
11220971522cSBarry Smith 
1123704ba839SBarry Smith EXTERN_C_BEGIN
1124704ba839SBarry Smith #undef __FUNCT__
1125704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
11267087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1127704ba839SBarry Smith {
1128704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1129704ba839SBarry Smith   PetscErrorCode    ierr;
1130704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1131704ba839SBarry Smith   char              prefix[128];
1132704ba839SBarry Smith 
1133704ba839SBarry Smith   PetscFunctionBegin;
11346c924f48SJed Brown   if (jac->splitdefined) {
11356c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
11366c924f48SJed Brown     PetscFunctionReturn(0);
11376c924f48SJed Brown   }
113816913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1139a04f6461SBarry Smith   if (splitname) {
1140db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1141a04f6461SBarry Smith   } else {
1142b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1143b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1144a04f6461SBarry Smith   }
11451ab39975SBarry Smith   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1146b5787286SJed Brown   ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
1147b5787286SJed Brown   ilink->is      = is;
1148b5787286SJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1149b5787286SJed Brown   ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1150b5787286SJed Brown   ilink->is_col  = is;
1151704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1152704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
1153*20252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1154704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
11559005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1156704ba839SBarry Smith 
1157a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1158704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1159704ba839SBarry Smith 
1160704ba839SBarry Smith   if (!next) {
1161704ba839SBarry Smith     jac->head       = ilink;
1162704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1163704ba839SBarry Smith   } else {
1164704ba839SBarry Smith     while (next->next) {
1165704ba839SBarry Smith       next = next->next;
1166704ba839SBarry Smith     }
1167704ba839SBarry Smith     next->next      = ilink;
1168704ba839SBarry Smith     ilink->previous = next;
1169704ba839SBarry Smith   }
1170704ba839SBarry Smith   jac->nsplits++;
1171704ba839SBarry Smith 
1172704ba839SBarry Smith   PetscFunctionReturn(0);
1173704ba839SBarry Smith }
1174704ba839SBarry Smith EXTERN_C_END
1175704ba839SBarry Smith 
11760971522cSBarry Smith #undef __FUNCT__
11770971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
11780971522cSBarry Smith /*@
11790971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
11800971522cSBarry Smith 
1181ad4df100SBarry Smith     Logically Collective on PC
11820971522cSBarry Smith 
11830971522cSBarry Smith     Input Parameters:
11840971522cSBarry Smith +   pc  - the preconditioner context
1185a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
11860971522cSBarry Smith .   n - the number of fields in this split
1187db4c96c1SJed Brown -   fields - the fields in this split
11880971522cSBarry Smith 
11890971522cSBarry Smith     Level: intermediate
11900971522cSBarry Smith 
1191d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1192d32f9abdSBarry Smith 
11937287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1194d32f9abdSBarry 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
1195d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1196d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1197d32f9abdSBarry Smith 
1198db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1199db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1200db4c96c1SJed Brown 
12015d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
12025d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
12035d4c12cdSJungho Lee      available when this routine is called.
12045d4c12cdSJungho Lee 
1205d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
12060971522cSBarry Smith 
12070971522cSBarry Smith @*/
12085d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
12090971522cSBarry Smith {
12104ac538c5SBarry Smith   PetscErrorCode ierr;
12110971522cSBarry Smith 
12120971522cSBarry Smith   PetscFunctionBegin;
12130700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1214db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1215db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1216db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
12175d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
12180971522cSBarry Smith   PetscFunctionReturn(0);
12190971522cSBarry Smith }
12200971522cSBarry Smith 
12210971522cSBarry Smith #undef __FUNCT__
1222704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1223704ba839SBarry Smith /*@
1224704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1225704ba839SBarry Smith 
1226ad4df100SBarry Smith     Logically Collective on PC
1227704ba839SBarry Smith 
1228704ba839SBarry Smith     Input Parameters:
1229704ba839SBarry Smith +   pc  - the preconditioner context
1230a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1231db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1232704ba839SBarry Smith 
1233d32f9abdSBarry Smith 
1234a6ffb8dbSJed Brown     Notes:
1235a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1236a6ffb8dbSJed Brown 
1237db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1238db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1239d32f9abdSBarry Smith 
1240704ba839SBarry Smith     Level: intermediate
1241704ba839SBarry Smith 
1242704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1243704ba839SBarry Smith 
1244704ba839SBarry Smith @*/
12457087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1246704ba839SBarry Smith {
12474ac538c5SBarry Smith   PetscErrorCode ierr;
1248704ba839SBarry Smith 
1249704ba839SBarry Smith   PetscFunctionBegin;
12500700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12517b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1252db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
12534ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1254704ba839SBarry Smith   PetscFunctionReturn(0);
1255704ba839SBarry Smith }
1256704ba839SBarry Smith 
1257704ba839SBarry Smith #undef __FUNCT__
125857a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
125957a9adfeSMatthew G Knepley /*@
126057a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
126157a9adfeSMatthew G Knepley 
126257a9adfeSMatthew G Knepley     Logically Collective on PC
126357a9adfeSMatthew G Knepley 
126457a9adfeSMatthew G Knepley     Input Parameters:
126557a9adfeSMatthew G Knepley +   pc  - the preconditioner context
126657a9adfeSMatthew G Knepley -   splitname - name of this split
126757a9adfeSMatthew G Knepley 
126857a9adfeSMatthew G Knepley     Output Parameter:
126957a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
127057a9adfeSMatthew G Knepley 
127157a9adfeSMatthew G Knepley     Level: intermediate
127257a9adfeSMatthew G Knepley 
127357a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
127457a9adfeSMatthew G Knepley 
127557a9adfeSMatthew G Knepley @*/
127657a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
127757a9adfeSMatthew G Knepley {
127857a9adfeSMatthew G Knepley   PetscErrorCode ierr;
127957a9adfeSMatthew G Knepley 
128057a9adfeSMatthew G Knepley   PetscFunctionBegin;
128157a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
128257a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
128357a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
128457a9adfeSMatthew G Knepley   {
128557a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
128657a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
128757a9adfeSMatthew G Knepley     PetscBool         found;
128857a9adfeSMatthew G Knepley 
128957a9adfeSMatthew G Knepley     *is = PETSC_NULL;
129057a9adfeSMatthew G Knepley     while(ilink) {
129157a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
129257a9adfeSMatthew G Knepley       if (found) {
129357a9adfeSMatthew G Knepley         *is = ilink->is;
129457a9adfeSMatthew G Knepley         break;
129557a9adfeSMatthew G Knepley       }
129657a9adfeSMatthew G Knepley       ilink = ilink->next;
129757a9adfeSMatthew G Knepley     }
129857a9adfeSMatthew G Knepley   }
129957a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
130057a9adfeSMatthew G Knepley }
130157a9adfeSMatthew G Knepley 
130257a9adfeSMatthew G Knepley #undef __FUNCT__
130351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
130451f519a2SBarry Smith /*@
130551f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
130651f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
130751f519a2SBarry Smith 
1308ad4df100SBarry Smith     Logically Collective on PC
130951f519a2SBarry Smith 
131051f519a2SBarry Smith     Input Parameters:
131151f519a2SBarry Smith +   pc  - the preconditioner context
131251f519a2SBarry Smith -   bs - the block size
131351f519a2SBarry Smith 
131451f519a2SBarry Smith     Level: intermediate
131551f519a2SBarry Smith 
131651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
131751f519a2SBarry Smith 
131851f519a2SBarry Smith @*/
13197087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
132051f519a2SBarry Smith {
13214ac538c5SBarry Smith   PetscErrorCode ierr;
132251f519a2SBarry Smith 
132351f519a2SBarry Smith   PetscFunctionBegin;
13240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1325c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
13264ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
132751f519a2SBarry Smith   PetscFunctionReturn(0);
132851f519a2SBarry Smith }
132951f519a2SBarry Smith 
133051f519a2SBarry Smith #undef __FUNCT__
133169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
13320971522cSBarry Smith /*@C
133369a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
13340971522cSBarry Smith 
133569a612a9SBarry Smith    Collective on KSP
13360971522cSBarry Smith 
13370971522cSBarry Smith    Input Parameter:
13380971522cSBarry Smith .  pc - the preconditioner context
13390971522cSBarry Smith 
13400971522cSBarry Smith    Output Parameters:
134113e0d083SBarry Smith +  n - the number of splits
134269a612a9SBarry Smith -  pc - the array of KSP contexts
13430971522cSBarry Smith 
13440971522cSBarry Smith    Note:
1345d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1346d32f9abdSBarry Smith    (not the KSP just the array that contains them).
13470971522cSBarry Smith 
134869a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
13490971522cSBarry Smith 
1350196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1351196cc216SBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the
1352196cc216SBarry Smith       KSP array must be.
1353196cc216SBarry Smith 
1354196cc216SBarry Smith 
13550971522cSBarry Smith    Level: advanced
13560971522cSBarry Smith 
13570971522cSBarry Smith .seealso: PCFIELDSPLIT
13580971522cSBarry Smith @*/
13597087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
13600971522cSBarry Smith {
13614ac538c5SBarry Smith   PetscErrorCode ierr;
13620971522cSBarry Smith 
13630971522cSBarry Smith   PetscFunctionBegin;
13640700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
136513e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
13664ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
13670971522cSBarry Smith   PetscFunctionReturn(0);
13680971522cSBarry Smith }
13690971522cSBarry Smith 
1370e69d4d44SBarry Smith #undef __FUNCT__
1371e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1372e69d4d44SBarry Smith /*@
1373e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1374a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1375e69d4d44SBarry Smith 
1376e69d4d44SBarry Smith     Collective on PC
1377e69d4d44SBarry Smith 
1378e69d4d44SBarry Smith     Input Parameters:
1379e69d4d44SBarry Smith +   pc  - the preconditioner context
1380fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1381084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1382084e4875SJed Brown 
1383e69d4d44SBarry Smith     Options Database:
1384084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1385e69d4d44SBarry Smith 
1386fd1303e9SJungho Lee     Notes:
1387fd1303e9SJungho Lee $    If ptype is
1388fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1389fd1303e9SJungho Lee $             to this function).
1390fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1391fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1392fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1393fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1394fd1303e9SJungho Lee $             preconditioner
1395fd1303e9SJungho Lee 
1396fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1397fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1398fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1399fd1303e9SJungho Lee 
1400fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1401fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1402fd1303e9SJungho Lee 
1403e69d4d44SBarry Smith     Level: intermediate
1404e69d4d44SBarry Smith 
1405fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1406e69d4d44SBarry Smith 
1407e69d4d44SBarry Smith @*/
14087087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1409e69d4d44SBarry Smith {
14104ac538c5SBarry Smith   PetscErrorCode ierr;
1411e69d4d44SBarry Smith 
1412e69d4d44SBarry Smith   PetscFunctionBegin;
14130700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14144ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1415e69d4d44SBarry Smith   PetscFunctionReturn(0);
1416e69d4d44SBarry Smith }
1417e69d4d44SBarry Smith 
1418e69d4d44SBarry Smith EXTERN_C_BEGIN
1419e69d4d44SBarry Smith #undef __FUNCT__
1420e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
14217087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1422e69d4d44SBarry Smith {
1423e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1424084e4875SJed Brown   PetscErrorCode  ierr;
1425e69d4d44SBarry Smith 
1426e69d4d44SBarry Smith   PetscFunctionBegin;
1427084e4875SJed Brown   jac->schurpre = ptype;
1428084e4875SJed Brown   if (pre) {
14296bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1430084e4875SJed Brown     jac->schur_user = pre;
1431084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1432084e4875SJed Brown   }
1433e69d4d44SBarry Smith   PetscFunctionReturn(0);
1434e69d4d44SBarry Smith }
1435e69d4d44SBarry Smith EXTERN_C_END
1436e69d4d44SBarry Smith 
143730ad9308SMatthew Knepley #undef __FUNCT__
1438c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1439ab1df9f5SJed Brown /*@
1440c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1441ab1df9f5SJed Brown 
1442ab1df9f5SJed Brown     Collective on PC
1443ab1df9f5SJed Brown 
1444ab1df9f5SJed Brown     Input Parameters:
1445ab1df9f5SJed Brown +   pc  - the preconditioner context
1446c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1447ab1df9f5SJed Brown 
1448ab1df9f5SJed Brown     Options Database:
1449c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1450ab1df9f5SJed Brown 
1451ab1df9f5SJed Brown 
1452ab1df9f5SJed Brown     Level: intermediate
1453ab1df9f5SJed Brown 
1454ab1df9f5SJed Brown     Notes:
1455ab1df9f5SJed Brown     The FULL factorization is
1456ab1df9f5SJed Brown 
1457ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1458ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1459ab1df9f5SJed Brown 
14606be4592eSBarry 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,
14616be4592eSBarry 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).
1462ab1df9f5SJed Brown 
14636be4592eSBarry 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
14646be4592eSBarry 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
14656be4592eSBarry 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,
14666be4592eSBarry 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.
1467ab1df9f5SJed Brown 
14686be4592eSBarry 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
14696be4592eSBarry 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).
1470ab1df9f5SJed Brown 
1471ab1df9f5SJed Brown     References:
1472ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1473ab1df9f5SJed Brown 
1474ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1475ab1df9f5SJed Brown 
1476ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1477ab1df9f5SJed Brown @*/
1478c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1479ab1df9f5SJed Brown {
1480ab1df9f5SJed Brown   PetscErrorCode ierr;
1481ab1df9f5SJed Brown 
1482ab1df9f5SJed Brown   PetscFunctionBegin;
1483ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1484c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1485ab1df9f5SJed Brown   PetscFunctionReturn(0);
1486ab1df9f5SJed Brown }
1487ab1df9f5SJed Brown 
1488ab1df9f5SJed Brown #undef __FUNCT__
1489c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1490c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1491ab1df9f5SJed Brown {
1492ab1df9f5SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1493ab1df9f5SJed Brown 
1494ab1df9f5SJed Brown   PetscFunctionBegin;
1495ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1496ab1df9f5SJed Brown   PetscFunctionReturn(0);
1497ab1df9f5SJed Brown }
1498ab1df9f5SJed Brown 
1499ab1df9f5SJed Brown #undef __FUNCT__
150030ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
150130ad9308SMatthew Knepley /*@C
15028c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
150330ad9308SMatthew Knepley 
150430ad9308SMatthew Knepley    Collective on KSP
150530ad9308SMatthew Knepley 
150630ad9308SMatthew Knepley    Input Parameter:
150730ad9308SMatthew Knepley .  pc - the preconditioner context
150830ad9308SMatthew Knepley 
150930ad9308SMatthew Knepley    Output Parameters:
1510a04f6461SBarry Smith +  A00 - the (0,0) block
1511a04f6461SBarry Smith .  A01 - the (0,1) block
1512a04f6461SBarry Smith .  A10 - the (1,0) block
1513a04f6461SBarry Smith -  A11 - the (1,1) block
151430ad9308SMatthew Knepley 
151530ad9308SMatthew Knepley    Level: advanced
151630ad9308SMatthew Knepley 
151730ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
151830ad9308SMatthew Knepley @*/
1519a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
152030ad9308SMatthew Knepley {
152130ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
152230ad9308SMatthew Knepley 
152330ad9308SMatthew Knepley   PetscFunctionBegin;
15240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1525c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1526a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1527a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1528a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1529a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
153030ad9308SMatthew Knepley   PetscFunctionReturn(0);
153130ad9308SMatthew Knepley }
153230ad9308SMatthew Knepley 
153379416396SBarry Smith EXTERN_C_BEGIN
153479416396SBarry Smith #undef __FUNCT__
153579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
15367087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
153779416396SBarry Smith {
153879416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1539e69d4d44SBarry Smith   PetscErrorCode ierr;
154079416396SBarry Smith 
154179416396SBarry Smith   PetscFunctionBegin;
154279416396SBarry Smith   jac->type = type;
15433b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
15443b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
15453b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1546e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1547e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1548c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1549e69d4d44SBarry Smith 
15503b224e63SBarry Smith   } else {
15513b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
15523b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1553e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
15549e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1555c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
15563b224e63SBarry Smith   }
155779416396SBarry Smith   PetscFunctionReturn(0);
155879416396SBarry Smith }
155979416396SBarry Smith EXTERN_C_END
156079416396SBarry Smith 
156151f519a2SBarry Smith EXTERN_C_BEGIN
156251f519a2SBarry Smith #undef __FUNCT__
156351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
15647087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
156551f519a2SBarry Smith {
156651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
156751f519a2SBarry Smith 
156851f519a2SBarry Smith   PetscFunctionBegin;
156965e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
157065e19b50SBarry 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);
157151f519a2SBarry Smith   jac->bs = bs;
157251f519a2SBarry Smith   PetscFunctionReturn(0);
157351f519a2SBarry Smith }
157451f519a2SBarry Smith EXTERN_C_END
157551f519a2SBarry Smith 
157679416396SBarry Smith #undef __FUNCT__
157779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1578bc08b0f1SBarry Smith /*@
157979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
158079416396SBarry Smith 
158179416396SBarry Smith    Collective on PC
158279416396SBarry Smith 
158379416396SBarry Smith    Input Parameter:
158479416396SBarry Smith .  pc - the preconditioner context
158581540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
158679416396SBarry Smith 
158779416396SBarry Smith    Options Database Key:
1588a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
158979416396SBarry Smith 
1590b02e2d75SMatthew G Knepley    Level: Intermediate
159179416396SBarry Smith 
159279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
159379416396SBarry Smith 
159479416396SBarry Smith .seealso: PCCompositeSetType()
159579416396SBarry Smith 
159679416396SBarry Smith @*/
15977087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
159879416396SBarry Smith {
15994ac538c5SBarry Smith   PetscErrorCode ierr;
160079416396SBarry Smith 
160179416396SBarry Smith   PetscFunctionBegin;
16020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16034ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
160479416396SBarry Smith   PetscFunctionReturn(0);
160579416396SBarry Smith }
160679416396SBarry Smith 
1607b02e2d75SMatthew G Knepley #undef __FUNCT__
1608b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType"
1609b02e2d75SMatthew G Knepley /*@
1610b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
1611b02e2d75SMatthew G Knepley 
1612b02e2d75SMatthew G Knepley   Not collective
1613b02e2d75SMatthew G Knepley 
1614b02e2d75SMatthew G Knepley   Input Parameter:
1615b02e2d75SMatthew G Knepley . pc - the preconditioner context
1616b02e2d75SMatthew G Knepley 
1617b02e2d75SMatthew G Knepley   Output Parameter:
1618b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1619b02e2d75SMatthew G Knepley 
1620b02e2d75SMatthew G Knepley   Level: Intermediate
1621b02e2d75SMatthew G Knepley 
1622b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1623b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
1624b02e2d75SMatthew G Knepley @*/
1625b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
1626b02e2d75SMatthew G Knepley {
1627b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
1628b02e2d75SMatthew G Knepley 
1629b02e2d75SMatthew G Knepley   PetscFunctionBegin;
1630b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1631b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
1632b02e2d75SMatthew G Knepley   *type = jac->type;
1633b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
1634b02e2d75SMatthew G Knepley }
1635b02e2d75SMatthew G Knepley 
16360971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
16370971522cSBarry Smith /*MC
1638a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1639a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
16400971522cSBarry Smith 
1641edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1642edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
16430971522cSBarry Smith 
1644a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
164569a612a9SBarry Smith          and set the options directly on the resulting KSP object
16460971522cSBarry Smith 
16470971522cSBarry Smith    Level: intermediate
16480971522cSBarry Smith 
164979416396SBarry Smith    Options Database Keys:
165081540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
165181540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
165281540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
165381540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
16540f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
16550f188ba9SJed Brown .   -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag
1656435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
165779416396SBarry Smith 
16585d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
16595d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
16605d4c12cdSJungho Lee 
1661c8a0d604SMatthew G Knepley    Notes:
1662c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1663d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1664d32f9abdSBarry Smith 
1665d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1666d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1667d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1668d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1669d32f9abdSBarry Smith 
1670c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1671c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1672c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1673c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1674c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1675a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1676c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1677c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1678c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1679a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1680c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1681c9c6ffaaSJed Brown      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
16825668aaf4SBarry Smith      diag gives
1683c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1684c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
16855668aaf4SBarry 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
1686c8a0d604SMatthew G Knepley $              (  A00   0 )
1687c8a0d604SMatthew G Knepley $              (  A10   S )
1688c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1689c8a0d604SMatthew G Knepley $              ( A00 A01 )
1690c8a0d604SMatthew G Knepley $              (  0   S  )
1691c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1692e69d4d44SBarry Smith 
1693edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1694edf189efSBarry Smith      is used automatically for a second block.
1695edf189efSBarry Smith 
1696ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1697ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1698ff218e97SBarry Smith 
1699ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1700ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1701ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
17020716a85fSBarry Smith 
1703a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
17040971522cSBarry Smith 
17057e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1706e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
17070971522cSBarry Smith M*/
17080971522cSBarry Smith 
17090971522cSBarry Smith EXTERN_C_BEGIN
17100971522cSBarry Smith #undef __FUNCT__
17110971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
17127087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
17130971522cSBarry Smith {
17140971522cSBarry Smith   PetscErrorCode ierr;
17150971522cSBarry Smith   PC_FieldSplit  *jac;
17160971522cSBarry Smith 
17170971522cSBarry Smith   PetscFunctionBegin;
171838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
17190971522cSBarry Smith   jac->bs        = -1;
17200971522cSBarry Smith   jac->nsplits   = 0;
17213e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1722e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1723c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
172451f519a2SBarry Smith 
17250971522cSBarry Smith   pc->data     = (void*)jac;
17260971522cSBarry Smith 
17270971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1728421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
17290971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1730574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
17310971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
17320971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
17330971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
17340971522cSBarry Smith   pc->ops->applyrichardson   = 0;
17350971522cSBarry Smith 
173669a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
173769a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
17380971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
17390971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1740704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1741704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
174279416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
174379416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
174451f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
174551f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
17460971522cSBarry Smith   PetscFunctionReturn(0);
17470971522cSBarry Smith }
17480971522cSBarry Smith EXTERN_C_END
17490971522cSBarry Smith 
17500971522cSBarry Smith 
1751a541d17aSBarry Smith 
1752