xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
1dba47a55SKris Buschelman 
2f5f0d762SBarry Smith 
3af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
4a80b646eSBarry Smith #include <petsc/private/kspimpl.h>    /*  This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
51e25c274SJed Brown #include <petscdm.h>
60971522cSBarry Smith 
72e71c61dSMatthew G. Knepley const char *const PCFieldSplitSchurPreTypes[] = {"SELF","SELFP","A11","USER","FULL","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
8c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
9c5d2311dSJed Brown 
109df09d43SBarry Smith PetscLogEvent KSP_Solve_FS_0,KSP_Solve_FS_1,KSP_Solve_FS_S,KSP_Solve_FS_U,KSP_Solve_FS_L,KSP_Solve_FS_2,KSP_Solve_FS_3,KSP_Solve_FS_4;
119df09d43SBarry Smith 
120971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
130971522cSBarry Smith struct _PC_FieldSplitLink {
1469a612a9SBarry Smith   KSP               ksp;
15443836d0SMatthew G Knepley   Vec               x,y,z;
16db4c96c1SJed Brown   char              *splitname;
170971522cSBarry Smith   PetscInt          nfields;
185d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
191b9fc7fcSBarry Smith   VecScatter        sctx;
20f5f0d762SBarry Smith   IS                is,is_col;
2151f519a2SBarry Smith   PC_FieldSplitLink next,previous;
229df09d43SBarry Smith   PetscLogEvent     event;
230971522cSBarry Smith };
240971522cSBarry Smith 
250971522cSBarry Smith typedef struct {
2668dd23aaSBarry Smith   PCCompositeType type;
27ace3abfcSBarry Smith   PetscBool       defaultsplit;                    /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
28ace3abfcSBarry Smith   PetscBool       splitdefined;                    /* Flag is set after the splits have been defined, to prevent more splits from being added */
2930ad9308SMatthew Knepley   PetscInt        bs;                              /* Block size for IS and Mat structures */
3030ad9308SMatthew Knepley   PetscInt        nsplits;                         /* Number of field divisions defined */
3179416396SBarry Smith   Vec             *x,*y,w1,w2;
32519d70e2SJed Brown   Mat             *mat;                            /* The diagonal block for each split */
33519d70e2SJed Brown   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
3430ad9308SMatthew Knepley   Mat             *Afield;                         /* The rows of the matrix associated with each split */
35ace3abfcSBarry Smith   PetscBool       issetup;
362fa5cd67SKarl Rupp 
3730ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3830ad9308SMatthew Knepley   Mat                       B;                     /* The (0,1) block */
3930ad9308SMatthew Knepley   Mat                       C;                     /* The (1,0) block */
40443836d0SMatthew G Knepley   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
41a7476a74SDmitry Karpeev   Mat                       schurp;                /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
42084e4875SJed Brown   Mat                       schur_user;            /* User-provided preconditioning matrix for the Schur complement */
43084e4875SJed Brown   PCFieldSplitSchurPreType  schurpre;              /* Determines which preconditioning matrix is used for the Schur complement */
44c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
4530ad9308SMatthew Knepley   KSP                       kspschur;              /* The solver for S */
46443836d0SMatthew G Knepley   KSP                       kspupper;              /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
47c096484dSStefano Zampini   PetscScalar               schurscale;            /* Scaling factor for the Schur complement solution with DIAG factorization */
48c096484dSStefano Zampini 
4997bbdb24SBarry Smith   PC_FieldSplitLink         head;
506dbb499eSCian Wilson   PetscBool                 isrestrict;             /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
51c1570756SJed Brown   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
524ab8060aSDmitry Karpeev   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
53c84da90fSDmitry Karpeev   PetscBool                 diag_use_amat;          /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
54c84da90fSDmitry Karpeev   PetscBool                 offdiag_use_amat;       /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
550971522cSBarry Smith } PC_FieldSplit;
560971522cSBarry Smith 
5716913363SBarry Smith /*
58*95452b02SPatrick Sanan     Notes:
59*95452b02SPatrick Sanan     there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
6016913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
6116913363SBarry Smith    PC you could change this.
6216913363SBarry Smith */
63084e4875SJed Brown 
64e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
65084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
66084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
67084e4875SJed Brown {
68084e4875SJed Brown   switch (jac->schurpre) {
69084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
70a7476a74SDmitry Karpeev   case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
71e87fab1bSBarry Smith   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
72e74569cdSMatthew G. Knepley   case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
73084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
74084e4875SJed Brown   default:
75084e4875SJed Brown     return jac->schur_user ? jac->schur_user : jac->pmat[1];
76084e4875SJed Brown   }
77084e4875SJed Brown }
78084e4875SJed Brown 
79084e4875SJed Brown 
809804daf3SBarry Smith #include <petscdraw.h>
810971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
820971522cSBarry Smith {
830971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
840971522cSBarry Smith   PetscErrorCode    ierr;
85d9884438SBarry Smith   PetscBool         iascii,isdraw;
860971522cSBarry Smith   PetscInt          i,j;
875a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
880971522cSBarry Smith 
890971522cSBarry Smith   PetscFunctionBegin;
90251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
91d9884438SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
920971522cSBarry Smith   if (iascii) {
932c9966d7SBarry Smith     if (jac->bs > 0) {
9451f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
952c9966d7SBarry Smith     } else {
962c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
972c9966d7SBarry Smith     }
98f5236f50SJed Brown     if (pc->useAmat) {
99f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
100a3df900dSMatthew G Knepley     }
101c84da90fSDmitry Karpeev     if (jac->diag_use_amat) {
102c84da90fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for diagonal blocks\n");CHKERRQ(ierr);
103c84da90fSDmitry Karpeev     }
104c84da90fSDmitry Karpeev     if (jac->offdiag_use_amat) {
105c84da90fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for off-diagonal blocks\n");CHKERRQ(ierr);
106c84da90fSDmitry Karpeev     }
10769a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
1080971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1090971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
1101ab39975SBarry Smith       if (ilink->fields) {
1110971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
11279416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1135a9f2f41SSatish Balay         for (j=0; j<ilink->nfields; j++) {
11479416396SBarry Smith           if (j > 0) {
11579416396SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
11679416396SBarry Smith           }
1175a9f2f41SSatish Balay           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1180971522cSBarry Smith         }
1190971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
12079416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1211ab39975SBarry Smith       } else {
1221ab39975SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1231ab39975SBarry Smith       }
1245a9f2f41SSatish Balay       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1255a9f2f41SSatish Balay       ilink = ilink->next;
1260971522cSBarry Smith     }
1270971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1282fa5cd67SKarl Rupp   }
1292fa5cd67SKarl Rupp 
1302fa5cd67SKarl Rupp  if (isdraw) {
131d9884438SBarry Smith     PetscDraw draw;
132d9884438SBarry Smith     PetscReal x,y,w,wd;
133d9884438SBarry Smith 
134d9884438SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
135d9884438SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
136d9884438SBarry Smith     w    = 2*PetscMin(1.0 - x,x);
137d9884438SBarry Smith     wd   = w/(jac->nsplits + 1);
138d9884438SBarry Smith     x    = x - wd*(jac->nsplits-1)/2.0;
139d9884438SBarry Smith     for (i=0; i<jac->nsplits; i++) {
140d9884438SBarry Smith       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
141d9884438SBarry Smith       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
142d9884438SBarry Smith       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
143d9884438SBarry Smith       x    += wd;
144d9884438SBarry Smith       ilink = ilink->next;
145d9884438SBarry Smith     }
1460971522cSBarry Smith   }
1470971522cSBarry Smith   PetscFunctionReturn(0);
1480971522cSBarry Smith }
1490971522cSBarry Smith 
1503b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1513b224e63SBarry Smith {
1523b224e63SBarry Smith   PC_FieldSplit              *jac = (PC_FieldSplit*)pc->data;
1533b224e63SBarry Smith   PetscErrorCode             ierr;
1544996c5bdSBarry Smith   PetscBool                  iascii,isdraw;
1553b224e63SBarry Smith   PetscInt                   i,j;
1563b224e63SBarry Smith   PC_FieldSplitLink          ilink = jac->head;
157a9908d51SBarry Smith   MatSchurComplementAinvType atype;
1583b224e63SBarry Smith 
1593b224e63SBarry Smith   PetscFunctionBegin;
160251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1614996c5bdSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1623b224e63SBarry Smith   if (iascii) {
1632c9966d7SBarry Smith     if (jac->bs > 0) {
164c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1652c9966d7SBarry Smith     } else {
1662c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1672c9966d7SBarry Smith     }
168f5236f50SJed Brown     if (pc->useAmat) {
169f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
170a3df900dSMatthew G Knepley     }
1713e8b8b31SMatthew G Knepley     switch (jac->schurpre) {
1723e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
173a9908d51SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);
174a9908d51SBarry Smith       break;
175a7476a74SDmitry Karpeev     case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
176a9908d51SBarry Smith       ierr = MatSchurComplementGetAinvType(jac->schur,&atype);CHKERRQ(ierr);
177d0a9c1a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses A00's %sdiagonal's inverse\n",atype == MAT_SCHUR_COMPLEMENT_AINV_DIAG ? "" : (atype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG ? "block " : "lumped "));CHKERRQ(ierr);break;
178e87fab1bSBarry Smith     case PC_FIELDSPLIT_SCHUR_PRE_A11:
179a9908d51SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
180a9908d51SBarry Smith       break;
181e74569cdSMatthew G. Knepley     case PC_FIELDSPLIT_SCHUR_PRE_FULL:
182a9908d51SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the exact Schur complement\n");CHKERRQ(ierr);
183a9908d51SBarry Smith       break;
1843e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1853e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1863e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1873e8b8b31SMatthew G Knepley       } else {
188e87fab1bSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
1893e8b8b31SMatthew G Knepley       }
1903e8b8b31SMatthew G Knepley       break;
1913e8b8b31SMatthew G Knepley     default:
19282f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1933e8b8b31SMatthew G Knepley     }
1943b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1953b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1963b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1973b224e63SBarry Smith       if (ilink->fields) {
1983b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1993b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
2003b224e63SBarry Smith         for (j=0; j<ilink->nfields; j++) {
2013b224e63SBarry Smith           if (j > 0) {
2023b224e63SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
2033b224e63SBarry Smith           }
2043b224e63SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
2053b224e63SBarry Smith         }
2063b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
2073b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
2083b224e63SBarry Smith       } else {
2093b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
2103b224e63SBarry Smith       }
2113b224e63SBarry Smith       ilink = ilink->next;
2123b224e63SBarry Smith     }
213435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");CHKERRQ(ierr);
2143b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
21506de4afeSJed Brown     if (jac->head) {
216443836d0SMatthew G Knepley       ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
21706de4afeSJed Brown     } else  {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
2183b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
21906de4afeSJed Brown     if (jac->head && jac->kspupper != jac->head->ksp) {
220443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr);
221443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
222b2750c55SJed Brown       if (jac->kspupper) {ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);}
223b2750c55SJed Brown       else {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
224443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
225443836d0SMatthew G Knepley     }
226435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
2273b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
22812cae6f2SJed Brown     if (jac->kspschur) {
2293b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
23012cae6f2SJed Brown     } else {
23112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
23212cae6f2SJed Brown     }
2333b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2343b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
23506de4afeSJed Brown   } else if (isdraw && jac->head) {
2364996c5bdSBarry Smith     PetscDraw draw;
2374996c5bdSBarry Smith     PetscReal x,y,w,wd,h;
2384996c5bdSBarry Smith     PetscInt  cnt = 2;
2394996c5bdSBarry Smith     char      str[32];
2404996c5bdSBarry Smith 
2414996c5bdSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2424996c5bdSBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
243c74581afSBarry Smith     if (jac->kspupper != jac->head->ksp) cnt++;
244c74581afSBarry Smith     w  = 2*PetscMin(1.0 - x,x);
245c74581afSBarry Smith     wd = w/(cnt + 1);
246c74581afSBarry Smith 
2474996c5bdSBarry Smith     ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
24851fa3d41SBarry Smith     ierr = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2494996c5bdSBarry Smith     y   -= h;
2504996c5bdSBarry Smith     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
251e87fab1bSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr);
2523b224e63SBarry Smith     } else {
2534996c5bdSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr);
2544996c5bdSBarry Smith     }
25551fa3d41SBarry Smith     ierr = PetscDrawStringBoxed(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2564996c5bdSBarry Smith     y   -= h;
2574996c5bdSBarry Smith     x    = x - wd*(cnt-1)/2.0;
2584996c5bdSBarry Smith 
2594996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2604996c5bdSBarry Smith     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
2614996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2624996c5bdSBarry Smith     if (jac->kspupper != jac->head->ksp) {
2634996c5bdSBarry Smith       x   += wd;
2644996c5bdSBarry Smith       ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2654996c5bdSBarry Smith       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
2664996c5bdSBarry Smith       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2674996c5bdSBarry Smith     }
2684996c5bdSBarry Smith     x   += wd;
2694996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2704996c5bdSBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
2714996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2723b224e63SBarry Smith   }
2733b224e63SBarry Smith   PetscFunctionReturn(0);
2743b224e63SBarry Smith }
2753b224e63SBarry Smith 
2766c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
2776c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
2786c924f48SJed Brown {
2796c924f48SJed Brown   PetscErrorCode ierr;
2806c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2815d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
2825d4c12cdSJungho Lee   PetscBool      flg,flg_col;
2835d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
2846c924f48SJed Brown 
2856c924f48SJed Brown   PetscFunctionBegin;
286785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields);CHKERRQ(ierr);
287785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields_col);CHKERRQ(ierr);
2886c924f48SJed Brown   for (i=0,flg=PETSC_TRUE;; i++) {
2898caf3d72SBarry Smith     ierr        = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
2908caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
2918caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
2926c924f48SJed Brown     nfields     = jac->bs;
29329499fbbSJungho Lee     nfields_col = jac->bs;
294c5929fdfSBarry Smith     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
295c5929fdfSBarry Smith     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
2966c924f48SJed Brown     if (!flg) break;
2975d4c12cdSJungho Lee     else if (flg && !flg_col) {
2985d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2995d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
3002fa5cd67SKarl Rupp     } else {
3015d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
3025d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
3035d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
3045d4c12cdSJungho Lee     }
3056c924f48SJed Brown   }
3066c924f48SJed Brown   if (i > 0) {
3076c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
3086c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
3096c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
3106c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
3116c924f48SJed Brown   }
3126c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
3135d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
3146c924f48SJed Brown   PetscFunctionReturn(0);
3156c924f48SJed Brown }
3166c924f48SJed Brown 
31769a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
3180971522cSBarry Smith {
3190971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3200971522cSBarry Smith   PetscErrorCode    ierr;
3215a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
3223a062f41SBarry Smith   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE,coupling = PETSC_FALSE;
3236c924f48SJed Brown   PetscInt          i;
3240971522cSBarry Smith 
3250971522cSBarry Smith   PetscFunctionBegin;
3267287d2a3SDmitry Karpeev   /*
327f5f0d762SBarry Smith    Kinda messy, but at least this now uses DMCreateFieldDecomposition().
3287287d2a3SDmitry Karpeev    Should probably be rewritten.
3297287d2a3SDmitry Karpeev    */
330f5f0d762SBarry Smith   if (!ilink) {
331c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
332c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_coupling",&coupling,NULL);CHKERRQ(ierr);
3333a062f41SBarry Smith     if (pc->dm && jac->dm_splits && !stokes && !coupling) {
334bafc1b83SMatthew G Knepley       PetscInt  numFields, f, i, j;
3350784a22cSJed Brown       char      **fieldNames;
3367b62db95SJungho Lee       IS        *fields;
337e7c4fc90SDmitry Karpeev       DM        *dms;
338bafc1b83SMatthew G Knepley       DM        subdm[128];
339bafc1b83SMatthew G Knepley       PetscBool flg;
340bafc1b83SMatthew G Knepley 
34116621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
342bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
343bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE;; i++) {
344bafc1b83SMatthew G Knepley         PetscInt ifields[128];
345bafc1b83SMatthew G Knepley         IS       compField;
346bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
347bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
348bafc1b83SMatthew G Knepley 
3498caf3d72SBarry Smith         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
350c5929fdfSBarry Smith         ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
351bafc1b83SMatthew G Knepley         if (!flg) break;
35282f516ccSBarry Smith         if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
353bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
354bafc1b83SMatthew G Knepley         if (nfields == 1) {
355bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
356bafc1b83SMatthew G Knepley         } else {
3578caf3d72SBarry Smith           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
358bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
3597287d2a3SDmitry Karpeev         }
360bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
361bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
362bafc1b83SMatthew G Knepley           f    = ifields[j];
3637b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
3647b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
3657b62db95SJungho Lee         }
366bafc1b83SMatthew G Knepley       }
367bafc1b83SMatthew G Knepley       if (i == 0) {
368bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
369bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
370bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
371bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
372bafc1b83SMatthew G Knepley         }
373bafc1b83SMatthew G Knepley       } else {
374d724dfffSBarry Smith         for (j=0; j<numFields; j++) {
375d724dfffSBarry Smith           ierr = DMDestroy(dms+j);CHKERRQ(ierr);
376d724dfffSBarry Smith         }
377d724dfffSBarry Smith         ierr = PetscFree(dms);CHKERRQ(ierr);
378785e854fSJed Brown         ierr = PetscMalloc1(i, &dms);CHKERRQ(ierr);
3792fa5cd67SKarl Rupp         for (j = 0; j < i; ++j) dms[j] = subdm[j];
380bafc1b83SMatthew G Knepley       }
3817b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
3827b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
383e7c4fc90SDmitry Karpeev       if (dms) {
3848b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
385bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
3867287d2a3SDmitry Karpeev           const char *prefix;
3877287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr);
3887287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr);
3897b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
3907b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
3917287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr);
392e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
3932fa5ba8aSJed Brown         }
3947b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
3958b8307b2SJed Brown       }
39666ffff09SJed Brown     } else {
397521d7252SBarry Smith       if (jac->bs <= 0) {
398704ba839SBarry Smith         if (pc->pmat) {
399521d7252SBarry Smith           ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
4002fa5cd67SKarl Rupp         } else jac->bs = 1;
401521d7252SBarry Smith       }
402d32f9abdSBarry Smith 
4036ce1633cSBarry Smith       if (stokes) {
4046ce1633cSBarry Smith         IS       zerodiags,rest;
4056ce1633cSBarry Smith         PetscInt nmin,nmax;
4066ce1633cSBarry Smith 
4076ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4086ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
4096ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
4106ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
4116ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
412fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
413fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4143a062f41SBarry Smith       } else if (coupling) {
4153a062f41SBarry Smith         IS       coupling,rest;
4163a062f41SBarry Smith         PetscInt nmin,nmax;
4173a062f41SBarry Smith 
4183a062f41SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4193a062f41SBarry Smith         ierr = MatFindOffBlockDiagonalEntries(pc->mat,&coupling);CHKERRQ(ierr);
4206bea0878SBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);CHKERRQ(ierr);
421e52d2c62SBarry Smith         ierr = ISSetIdentity(rest);CHKERRQ(ierr);
422d020c841SBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
423d020c841SBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",coupling);CHKERRQ(ierr);
4243a062f41SBarry Smith         ierr = ISDestroy(&coupling);CHKERRQ(ierr);
4253a062f41SBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4266ce1633cSBarry Smith       } else {
427c5929fdfSBarry Smith         ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
4287287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
429d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
430d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
4316c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
4326c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
433d32f9abdSBarry Smith         }
4346dbb499eSCian Wilson         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
435d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
436db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
4376c924f48SJed Brown             char splitname[8];
4388caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
4395d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
44079416396SBarry Smith           }
4415d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
442521d7252SBarry Smith         }
44366ffff09SJed Brown       }
4446ce1633cSBarry Smith     }
445edf189efSBarry Smith   } else if (jac->nsplits == 1) {
446edf189efSBarry Smith     if (ilink->is) {
447edf189efSBarry Smith       IS       is2;
448edf189efSBarry Smith       PetscInt nmin,nmax;
449edf189efSBarry Smith 
450edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
451edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
452db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
453fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
45482f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
455edf189efSBarry Smith   }
456d0af7cd3SBarry Smith 
45782f516ccSBarry Smith   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
45869a612a9SBarry Smith   PetscFunctionReturn(0);
45969a612a9SBarry Smith }
46069a612a9SBarry Smith 
461c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions,const char pre[], const char name[], char *value[], PetscBool *flg);
462514bf10dSMatthew G Knepley 
46369a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
46469a612a9SBarry Smith {
46569a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
46669a612a9SBarry Smith   PetscErrorCode    ierr;
4675a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
4682c9966d7SBarry Smith   PetscInt          i,nsplit;
4695f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
47069a612a9SBarry Smith 
47169a612a9SBarry Smith   PetscFunctionBegin;
47269a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
47397bbdb24SBarry Smith   nsplit = jac->nsplits;
4745a9f2f41SSatish Balay   ilink  = jac->head;
47597bbdb24SBarry Smith 
47697bbdb24SBarry Smith   /* get the matrices for each split */
477704ba839SBarry Smith   if (!jac->issetup) {
4781b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
47997bbdb24SBarry Smith 
480704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
481704ba839SBarry Smith 
4825d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4832c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4842c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
4852c9966d7SBarry Smith     }
48651f519a2SBarry Smith     bs     = jac->bs;
48797bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
4881b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4891b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4901b9fc7fcSBarry Smith       if (jac->defaultsplit) {
491ce94432eSBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
4925f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
493704ba839SBarry Smith       } else if (!ilink->is) {
494ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4955f4ab4e1SJungho Lee           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
496785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&ii);CHKERRQ(ierr);
497785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&jj);CHKERRQ(ierr);
4981b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
4991b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
5001b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
5015f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
50297bbdb24SBarry Smith             }
50397bbdb24SBarry Smith           }
504ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
505ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
50690e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is, nfields);CHKERRQ(ierr);
50790e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is_col, nfields);CHKERRQ(ierr);
508ccb205f8SBarry Smith         } else {
509ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
510ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
511ccb205f8SBarry Smith        }
5123e197d65SBarry Smith       }
513704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
5145f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
5155f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
516704ba839SBarry Smith       ilink = ilink->next;
5171b9fc7fcSBarry Smith     }
5181b9fc7fcSBarry Smith   }
5191b9fc7fcSBarry Smith 
520704ba839SBarry Smith   ilink = jac->head;
52197bbdb24SBarry Smith   if (!jac->pmat) {
522bdddcaaaSMatthew G Knepley     Vec xtmp;
523bdddcaaaSMatthew G Knepley 
5242a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
525785e854fSJed Brown     ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr);
526dcca6d9dSJed Brown     ierr = PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);CHKERRQ(ierr);
527cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
528bdddcaaaSMatthew G Knepley       MatNullSpace sp;
529bdddcaaaSMatthew G Knepley 
530a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
531a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
532a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
533a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
534a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
535a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
5362fa5cd67SKarl Rupp 
537a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
538a3df900dSMatthew G Knepley         }
539a3df900dSMatthew G Knepley       } else {
5403a062f41SBarry Smith         const char *prefix;
5417dae84e0SHong Zhang         ierr = MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
5423a062f41SBarry Smith         ierr = KSPGetOptionsPrefix(ilink->ksp,&prefix);CHKERRQ(ierr);
5433a062f41SBarry Smith         ierr = MatSetOptionsPrefix(jac->pmat[i],prefix);CHKERRQ(ierr);
5443a062f41SBarry Smith         ierr = MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");CHKERRQ(ierr);
545a3df900dSMatthew G Knepley       }
546bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
5472a7a6963SBarry Smith       ierr     = MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
5480298fd71SBarry Smith       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
549bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
5500298fd71SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
551ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
552ed1f0337SMatthew G Knepley       if (sp) {
553ed1f0337SMatthew G Knepley         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
554ed1f0337SMatthew G Knepley       }
555704ba839SBarry Smith       ilink = ilink->next;
556cf502942SBarry Smith     }
557bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
55897bbdb24SBarry Smith   } else {
559cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
560a3df900dSMatthew G Knepley       Mat pmat;
561a3df900dSMatthew G Knepley 
562a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
563a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
564a3df900dSMatthew G Knepley       if (!pmat) {
5657dae84e0SHong Zhang         ierr = MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
566a3df900dSMatthew G Knepley       }
567704ba839SBarry Smith       ilink = ilink->next;
568cf502942SBarry Smith     }
56997bbdb24SBarry Smith   }
5704e39094bSDmitry Karpeev   if (jac->diag_use_amat) {
571519d70e2SJed Brown     ilink = jac->head;
572519d70e2SJed Brown     if (!jac->mat) {
573785e854fSJed Brown       ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr);
574519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5757dae84e0SHong Zhang         ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
576519d70e2SJed Brown         ilink = ilink->next;
577519d70e2SJed Brown       }
578519d70e2SJed Brown     } else {
579519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5807dae84e0SHong Zhang         if (jac->mat[i]) {ierr = MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
581519d70e2SJed Brown         ilink = ilink->next;
582519d70e2SJed Brown       }
583519d70e2SJed Brown     }
584519d70e2SJed Brown   } else {
585519d70e2SJed Brown     jac->mat = jac->pmat;
586519d70e2SJed Brown   }
58797bbdb24SBarry Smith 
58853935eafSBarry Smith   /* Check for null space attached to IS */
58953935eafSBarry Smith   ilink = jac->head;
59053935eafSBarry Smith   for (i=0; i<nsplit; i++) {
59153935eafSBarry Smith     MatNullSpace sp;
59253935eafSBarry Smith 
59353935eafSBarry Smith     ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
59453935eafSBarry Smith     if (sp) {
59553935eafSBarry Smith       ierr = MatSetNullSpace(jac->mat[i], sp);CHKERRQ(ierr);
59653935eafSBarry Smith     }
59753935eafSBarry Smith     ilink = ilink->next;
59853935eafSBarry Smith   }
59953935eafSBarry Smith 
6006c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
60168dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
6024e39094bSDmitry Karpeev     /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
60368dd23aaSBarry Smith     ilink = jac->head;
604e52d2c62SBarry Smith     if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
605e52d2c62SBarry Smith       /* special case need where Afield[0] is not needed and only certain columns of Afield[1] are needed since update is only on those rows of the solution */
606e52d2c62SBarry Smith       if (!jac->Afield) {
607e52d2c62SBarry Smith         ierr = PetscCalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
60880c96bb1SFande Kong         if (jac->offdiag_use_amat) {
6097dae84e0SHong Zhang           ierr  = MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
610e52d2c62SBarry Smith         } else {
6117dae84e0SHong Zhang           ierr  = MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
61280c96bb1SFande Kong         }
61380c96bb1SFande Kong       } else {
61480c96bb1SFande Kong         if (jac->offdiag_use_amat) {
6157dae84e0SHong Zhang           ierr  = MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_REUSE_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
61680c96bb1SFande Kong         } else {
6177dae84e0SHong Zhang           ierr  = MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,MAT_REUSE_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
61880c96bb1SFande Kong         }
619e52d2c62SBarry Smith       }
620e52d2c62SBarry Smith     } else {
62168dd23aaSBarry Smith       if (!jac->Afield) {
622785e854fSJed Brown         ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
62368dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
62480c96bb1SFande Kong           if (jac->offdiag_use_amat) {
6257dae84e0SHong Zhang             ierr  = MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
62680c96bb1SFande Kong           } else {
6277dae84e0SHong Zhang             ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
62880c96bb1SFande Kong           }
62968dd23aaSBarry Smith           ilink = ilink->next;
63068dd23aaSBarry Smith         }
63168dd23aaSBarry Smith       } else {
63268dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
63380c96bb1SFande Kong           if (jac->offdiag_use_amat) {
6347dae84e0SHong Zhang             ierr  = MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63580c96bb1SFande Kong           } else {
6367dae84e0SHong Zhang             ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63780c96bb1SFande Kong           }
63868dd23aaSBarry Smith           ilink = ilink->next;
63968dd23aaSBarry Smith         }
64068dd23aaSBarry Smith       }
64168dd23aaSBarry Smith     }
642e52d2c62SBarry Smith   }
64368dd23aaSBarry Smith 
6443b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
6453b224e63SBarry Smith     IS          ccis;
646c096484dSStefano Zampini     PetscBool   isspd;
6474aa3045dSJed Brown     PetscInt    rstart,rend;
648093c86ffSJed Brown     char        lscname[256];
649093c86ffSJed Brown     PetscObject LSC_L;
650ce94432eSBarry Smith 
651ce94432eSBarry Smith     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
65268dd23aaSBarry Smith 
653c096484dSStefano Zampini     /* If pc->mat is SPD, don't scale by -1 the Schur complement */
654c096484dSStefano Zampini     if (jac->schurscale == (PetscScalar)-1.0) {
655c096484dSStefano Zampini       ierr = MatGetOption(pc->pmat,MAT_SPD,&isspd);CHKERRQ(ierr);
656c096484dSStefano Zampini       jac->schurscale = (isspd == PETSC_TRUE) ? 1.0 : -1.0;
657c096484dSStefano Zampini     }
658c096484dSStefano Zampini 
659e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
660e6cab6aaSJed Brown     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
661e6cab6aaSJed Brown 
6623b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
6633b224e63SBarry Smith     if (jac->schur) {
6640298fd71SBarry Smith       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
665443836d0SMatthew G Knepley 
666fb3147dbSMatthew G Knepley       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
6673b224e63SBarry Smith       ilink = jac->head;
66849bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6694e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6707dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
671c84da90fSDmitry Karpeev       } else {
6727dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
673c84da90fSDmitry Karpeev       }
674fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6753b224e63SBarry Smith       ilink = ilink->next;
67649bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6774e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6787dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
679c84da90fSDmitry Karpeev       } else {
6807dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
681c84da90fSDmitry Karpeev       }
682fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
683aa6c7ce3SBarry Smith       ierr  = MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
684a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
685a7476a74SDmitry Karpeev 	ierr = MatDestroy(&jac->schurp);CHKERRQ(ierr);
686a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
687a7476a74SDmitry Karpeev       }
688443836d0SMatthew G Knepley       if (kspA != kspInner) {
68923ee1639SBarry Smith         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
690443836d0SMatthew G Knepley       }
691443836d0SMatthew G Knepley       if (kspUpper != kspA) {
69223ee1639SBarry Smith         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
693443836d0SMatthew G Knepley       }
69423ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
6953b224e63SBarry Smith     } else {
696bafc1b83SMatthew G Knepley       const char   *Dprefix;
697470b340bSDmitry Karpeev       char         schurprefix[256], schurmatprefix[256];
698514bf10dSMatthew G Knepley       char         schurtestoption[256];
699bdddcaaaSMatthew G Knepley       MatNullSpace sp;
700514bf10dSMatthew G Knepley       PetscBool    flg;
7013b224e63SBarry Smith 
702a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
7033b224e63SBarry Smith       ilink = jac->head;
70449bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
7054e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
7067dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
707c84da90fSDmitry Karpeev       } else {
7087dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
709c84da90fSDmitry Karpeev       }
710fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
7113b224e63SBarry Smith       ilink = ilink->next;
71249bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
7134e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
7147dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
715c84da90fSDmitry Karpeev       } else {
7167dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
717c84da90fSDmitry Karpeev       }
718fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
71920252d06SBarry Smith 
720f5236f50SJed Brown       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
72120252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
72220252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
723bee83525SDmitry Karpeev       ierr = MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
724470b340bSDmitry Karpeev       ierr = PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
725470b340bSDmitry Karpeev       /* Note that the inner KSP is NOT going to inherit this prefix, and if it did, it would be reset just below.  Is that what we want? */
726470b340bSDmitry Karpeev       ierr = MatSetOptionsPrefix(jac->schur,schurmatprefix);CHKERRQ(ierr);
727470b340bSDmitry Karpeev       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
728895c21f2SBarry Smith       ierr = MatGetNullSpace(jac->mat[1], &sp);CHKERRQ(ierr);
72920252d06SBarry Smith       if (sp) {
73020252d06SBarry Smith         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
73120252d06SBarry Smith       }
73220252d06SBarry Smith 
73320252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
734c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
735514bf10dSMatthew G Knepley       if (flg) {
736514bf10dSMatthew G Knepley         DM  dmInner;
73721635b76SJed Brown         KSP kspInner;
73821635b76SJed Brown 
73921635b76SJed Brown         ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
74021635b76SJed Brown         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
74121635b76SJed Brown         /* Indent this deeper to emphasize the "inner" nature of this solver. */
74221635b76SJed Brown         ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr);
74321635b76SJed Brown         ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr);
744514bf10dSMatthew G Knepley 
745514bf10dSMatthew G Knepley         /* Set DM for new solver */
746bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
74721635b76SJed Brown         ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr);
74821635b76SJed Brown         ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr);
749514bf10dSMatthew G Knepley       } else {
75021635b76SJed Brown          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
75121635b76SJed Brown           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
75221635b76SJed Brown           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
75321635b76SJed Brown           * S = D - C A_inner^{-1} B, we expect S to be defined using an accurate definition of A_inner^{-1}, so we make
75421635b76SJed Brown           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
75521635b76SJed Brown           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
75621635b76SJed Brown         ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr);
757514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr);
758bafc1b83SMatthew G Knepley       }
75923ee1639SBarry Smith       ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
7605a9f2f41SSatish Balay       ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
76197bbdb24SBarry Smith       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
76297bbdb24SBarry Smith 
763443836d0SMatthew G Knepley       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
764c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
765443836d0SMatthew G Knepley       if (flg) {
766443836d0SMatthew G Knepley         DM dmInner;
767443836d0SMatthew G Knepley 
768443836d0SMatthew G Knepley         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
76982f516ccSBarry Smith         ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr);
770422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);CHKERRQ(ierr);
771443836d0SMatthew G Knepley         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
772443836d0SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
773443836d0SMatthew G Knepley         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
774443836d0SMatthew G Knepley         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
775443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
77623ee1639SBarry Smith         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
777443836d0SMatthew G Knepley         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
778443836d0SMatthew G Knepley       } else {
779443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
780443836d0SMatthew G Knepley         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
781443836d0SMatthew G Knepley       }
782443836d0SMatthew G Knepley 
783a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
784a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
785a7476a74SDmitry Karpeev       }
786ce94432eSBarry Smith       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
787422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);CHKERRQ(ierr);
78897bbdb24SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
78920252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
79097bbdb24SBarry Smith       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
7917233a360SDmitry Karpeev         PC pcschur;
7927233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
7937233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
79497bbdb24SBarry Smith         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
795e74569cdSMatthew G. Knepley       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
796e74569cdSMatthew G. Knepley         ierr = MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);CHKERRQ(ierr);
79797bbdb24SBarry Smith       }
79823ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
79997bbdb24SBarry Smith       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
80097bbdb24SBarry Smith       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
801c096484dSStefano Zampini       /* propagate DM */
802b20b4189SMatthew G. Knepley       {
803b20b4189SMatthew G. Knepley         DM sdm;
804b20b4189SMatthew G. Knepley         ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr);
805b20b4189SMatthew G. Knepley         if (sdm) {
806b20b4189SMatthew G. Knepley           ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr);
807b20b4189SMatthew G. Knepley           ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr);
808b20b4189SMatthew G. Knepley         }
809b20b4189SMatthew G. Knepley       }
81097bbdb24SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
81197bbdb24SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
81297bbdb24SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
81397bbdb24SBarry Smith     }
81497bbdb24SBarry Smith 
8155a9f2f41SSatish Balay     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
8168caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
817519d70e2SJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
8183b224e63SBarry Smith     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
819c1570756SJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
8208caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
82197bbdb24SBarry Smith     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
8225a9f2f41SSatish Balay     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
8230971522cSBarry Smith     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
82497bbdb24SBarry Smith   } else {
82568bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
8260971522cSBarry Smith     i     = 0;
8270971522cSBarry Smith     ilink = jac->head;
8280971522cSBarry Smith     while (ilink) {
82923ee1639SBarry Smith       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);CHKERRQ(ierr);
8300971522cSBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
8310971522cSBarry Smith       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
8320971522cSBarry Smith       i++;
8330971522cSBarry Smith       ilink = ilink->next;
8340971522cSBarry Smith     }
8353b224e63SBarry Smith   }
8363b224e63SBarry Smith 
837c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
8380971522cSBarry Smith   PetscFunctionReturn(0);
8390971522cSBarry Smith }
8400971522cSBarry Smith 
8415a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
842ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
843ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
8449df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
8455a9f2f41SSatish Balay    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
8469df09d43SBarry Smith    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
847ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
848ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
84979416396SBarry Smith 
8503b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
8513b224e63SBarry Smith {
8523b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8533b224e63SBarry Smith   PetscErrorCode    ierr;
8543b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
855443836d0SMatthew G Knepley   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
8563b224e63SBarry Smith 
8573b224e63SBarry Smith   PetscFunctionBegin;
858c5d2311dSJed Brown   switch (jac->schurfactorization) {
859c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
860a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
861c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
862c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
863c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8649df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
865443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8669df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
867c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
868c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8699df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
870c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8719df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
872c096484dSStefano Zampini     ierr = VecScale(ilinkD->y,jac->schurscale);CHKERRQ(ierr);
873c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
874c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
875c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
876c5d2311dSJed Brown     break;
877c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
878a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
879c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
880c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8819df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
882443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8839df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
884c5d2311dSJed Brown     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
885c5d2311dSJed Brown     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
886c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
887c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
888c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8899df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
890c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8919df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
892c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
893c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
894c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
895c5d2311dSJed Brown     break;
896c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
897a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
898c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
899c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9009df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
901c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
9029df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);    ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
903c5d2311dSJed Brown     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
904c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
905c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
906c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9079df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
908443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9099df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
910c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
911c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
912c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
913c5d2311dSJed Brown     break;
914c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
915a04f6461SBarry 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 */
9163b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9173b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9189df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
919443836d0SMatthew G Knepley     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9209df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
9213b224e63SBarry Smith     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
9223b224e63SBarry Smith     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
9233b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9243b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9253b224e63SBarry Smith 
9269df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
9273b224e63SBarry Smith     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
9289df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
9293b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9303b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9313b224e63SBarry Smith 
932443836d0SMatthew G Knepley     if (kspUpper == kspA) {
9333b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
9343b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
9359df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
936443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9379df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
938443836d0SMatthew G Knepley     } else {
9399df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
940443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9419df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
942443836d0SMatthew G Knepley       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
9439df09d43SBarry Smith       ierr = PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
944443836d0SMatthew G Knepley       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
9459df09d43SBarry Smith       ierr = PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
946443836d0SMatthew G Knepley       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
947443836d0SMatthew G Knepley     }
9483b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9493b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950c5d2311dSJed Brown   }
9513b224e63SBarry Smith   PetscFunctionReturn(0);
9523b224e63SBarry Smith }
9533b224e63SBarry Smith 
9540971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
9550971522cSBarry Smith {
9560971522cSBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
9570971522cSBarry Smith   PetscErrorCode     ierr;
9585a9f2f41SSatish Balay   PC_FieldSplitLink  ilink = jac->head;
959939b8a20SBarry Smith   PetscInt           cnt,bs;
960568ed31eSHong Zhang   KSPConvergedReason reason;
9610971522cSBarry Smith 
9620971522cSBarry Smith   PetscFunctionBegin;
96379416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
9641b9fc7fcSBarry Smith     if (jac->defaultsplit) {
965939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
966ce94432eSBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
967939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
968ce94432eSBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
9690971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
9705a9f2f41SSatish Balay       while (ilink) {
9719df09d43SBarry Smith         ierr  = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
9725a9f2f41SSatish Balay         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9739df09d43SBarry Smith         ierr  = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
974568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
975568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
976568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
977568ed31eSHong Zhang         }
9785a9f2f41SSatish Balay         ilink = ilink->next;
9790971522cSBarry Smith       }
9800971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
9811b9fc7fcSBarry Smith     } else {
982efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
9835a9f2f41SSatish Balay       while (ilink) {
9845a9f2f41SSatish Balay         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
985568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
986568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
987568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
988568ed31eSHong Zhang         }
9895a9f2f41SSatish Balay         ilink = ilink->next;
9901b9fc7fcSBarry Smith       }
9911b9fc7fcSBarry Smith     }
992e52d2c62SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
993e52d2c62SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
994e52d2c62SBarry Smith     /* solve on first block for first block variables */
995e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
996e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9979df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
998e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9999df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1000568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1001568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1002568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1003568ed31eSHong Zhang     }
1004e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1005e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1006e52d2c62SBarry Smith 
1007e52d2c62SBarry Smith     /* compute the residual only onto second block variables using first block variables */
1008e52d2c62SBarry Smith     ierr = MatMult(jac->Afield[1],ilink->y,ilink->next->x);CHKERRQ(ierr);
1009e52d2c62SBarry Smith     ilink = ilink->next;
1010e52d2c62SBarry Smith     ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
1011e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1012e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1013e52d2c62SBarry Smith 
1014e52d2c62SBarry Smith     /* solve on second block variables */
10159df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1016e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10179df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1018568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1019568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1020568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1021568ed31eSHong Zhang     }
1022e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1023e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
102416913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
102579416396SBarry Smith     if (!jac->w1) {
102679416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
102779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
102879416396SBarry Smith     }
1029efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
10305a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
1031568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1032568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1033568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1034568ed31eSHong Zhang     }
10353e197d65SBarry Smith     cnt  = 1;
10365a9f2f41SSatish Balay     while (ilink->next) {
10375a9f2f41SSatish Balay       ilink = ilink->next;
10383e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
10393e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
10403e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
10413e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10423e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10439df09d43SBarry Smith       ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
10443e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10459df09d43SBarry Smith       ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
10462618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
10472618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
10482618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
10492618eb8fSHong Zhang       }
10503e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10513e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10523e197d65SBarry Smith     }
105351f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
105411755939SBarry Smith       cnt -= 2;
105551f519a2SBarry Smith       while (ilink->previous) {
105651f519a2SBarry Smith         ilink = ilink->previous;
105711755939SBarry Smith         /* compute the residual only over the part of the vector needed */
105811755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
105911755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
106011755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
106111755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10629df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
106311755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10649df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1065568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1066568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1067568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
1068568ed31eSHong Zhang         }
106911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
107011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
107151f519a2SBarry Smith       }
107211755939SBarry Smith     }
1073ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
10740971522cSBarry Smith   PetscFunctionReturn(0);
10750971522cSBarry Smith }
10760971522cSBarry Smith 
1077421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1078ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1079ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
10809df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1081421e10b8SBarry Smith    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) ||                  \
10829df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1083ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1084ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
1085421e10b8SBarry Smith 
1086421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1087421e10b8SBarry Smith {
1088421e10b8SBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1089421e10b8SBarry Smith   PetscErrorCode     ierr;
1090421e10b8SBarry Smith   PC_FieldSplitLink  ilink = jac->head;
1091939b8a20SBarry Smith   PetscInt           bs;
1092291dd214SHong Zhang   KSPConvergedReason reason;
1093421e10b8SBarry Smith 
1094421e10b8SBarry Smith   PetscFunctionBegin;
1095421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1096421e10b8SBarry Smith     if (jac->defaultsplit) {
1097939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
1098ce94432eSBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1099939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
1100ce94432eSBarry Smith       if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1101421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
1102421e10b8SBarry Smith       while (ilink) {
11039df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1104421e10b8SBarry Smith         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
11059df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
11062618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11072618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11082618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
11092618eb8fSHong Zhang         }
1110421e10b8SBarry Smith         ilink = ilink->next;
1111421e10b8SBarry Smith       }
1112421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
1113421e10b8SBarry Smith     } else {
1114421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
1115421e10b8SBarry Smith       while (ilink) {
1116421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11172618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11182618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11192618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
11202618eb8fSHong Zhang         }
1121421e10b8SBarry Smith         ilink = ilink->next;
1122421e10b8SBarry Smith       }
1123421e10b8SBarry Smith     }
1124421e10b8SBarry Smith   } else {
1125421e10b8SBarry Smith     if (!jac->w1) {
1126421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
1127421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
1128421e10b8SBarry Smith     }
1129421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
1130421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1131421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11322618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11332618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11342618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
11352618eb8fSHong Zhang       }
1136421e10b8SBarry Smith       while (ilink->next) {
1137421e10b8SBarry Smith         ilink = ilink->next;
11389989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1139421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1140421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1141421e10b8SBarry Smith       }
1142421e10b8SBarry Smith       while (ilink->previous) {
1143421e10b8SBarry Smith         ilink = ilink->previous;
11449989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1145421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1146421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1147421e10b8SBarry Smith       }
1148421e10b8SBarry Smith     } else {
1149421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
1150421e10b8SBarry Smith         ilink = ilink->next;
1151421e10b8SBarry Smith       }
1152421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11532618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11542618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11552618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
11562618eb8fSHong Zhang       }
1157421e10b8SBarry Smith       while (ilink->previous) {
1158421e10b8SBarry Smith         ilink = ilink->previous;
11599989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1160421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1161421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1162421e10b8SBarry Smith       }
1163421e10b8SBarry Smith     }
1164421e10b8SBarry Smith   }
1165421e10b8SBarry Smith   PetscFunctionReturn(0);
1166421e10b8SBarry Smith }
1167421e10b8SBarry Smith 
1168574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
11690971522cSBarry Smith {
11700971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11710971522cSBarry Smith   PetscErrorCode    ierr;
11725a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
11730971522cSBarry Smith 
11740971522cSBarry Smith   PetscFunctionBegin;
11755a9f2f41SSatish Balay   while (ilink) {
1176f5f0d762SBarry Smith     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1177fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
1178fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
1179443836d0SMatthew G Knepley     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
1180fcfd50ebSBarry Smith     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
1181fcfd50ebSBarry Smith     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
1182b5787286SJed Brown     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1183f5f0d762SBarry Smith     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1184f5f0d762SBarry Smith     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
1185f5f0d762SBarry Smith     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
11865a9f2f41SSatish Balay     next  = ilink->next;
1187f5f0d762SBarry Smith     ierr  = PetscFree(ilink);CHKERRQ(ierr);
11885a9f2f41SSatish Balay     ilink = next;
11890971522cSBarry Smith   }
1190f5f0d762SBarry Smith   jac->head = NULL;
119105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1192574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
1193574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
1194574deadeSBarry Smith   } else if (jac->mat) {
11950298fd71SBarry Smith     jac->mat = NULL;
1196574deadeSBarry Smith   }
119797bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
119868dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
1199f5f0d762SBarry Smith   jac->nsplits = 0;
12006bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
12016bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
12026bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
1203a7476a74SDmitry Karpeev   ierr       = MatDestroy(&jac->schurp);CHKERRQ(ierr);
12046bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
12056bf464f9SBarry Smith   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1206d78dad28SBarry Smith   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
12076bf464f9SBarry Smith   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
12086bf464f9SBarry Smith   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
12096dbb499eSCian Wilson   jac->isrestrict = PETSC_FALSE;
1210574deadeSBarry Smith   PetscFunctionReturn(0);
1211574deadeSBarry Smith }
1212574deadeSBarry Smith 
1213574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1214574deadeSBarry Smith {
1215574deadeSBarry Smith   PetscErrorCode    ierr;
1216574deadeSBarry Smith 
1217574deadeSBarry Smith   PetscFunctionBegin;
1218574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1219c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1220bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr);
1221bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr);
1222bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr);
1223bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr);
1224bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr);
122529f8a81cSJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);CHKERRQ(ierr);
122637a82bf0SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);CHKERRQ(ierr);
1227bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr);
12286dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);CHKERRQ(ierr);
12290971522cSBarry Smith   PetscFunctionReturn(0);
12300971522cSBarry Smith }
12310971522cSBarry Smith 
12324416b707SBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
12330971522cSBarry Smith {
12341b9fc7fcSBarry Smith   PetscErrorCode  ierr;
12356c924f48SJed Brown   PetscInt        bs;
1236bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
12379dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
12383b224e63SBarry Smith   PCCompositeType ctype;
12391b9fc7fcSBarry Smith 
12400971522cSBarry Smith   PetscFunctionBegin;
1241e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"FieldSplit options");CHKERRQ(ierr);
12424ab8060aSDmitry Karpeev   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
124351f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
124451f519a2SBarry Smith   if (flg) {
124551f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
124651f519a2SBarry Smith   }
12472686e3e9SMatthew G. Knepley   jac->diag_use_amat = pc->useAmat;
12482686e3e9SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_fieldsplit_diag_use_amat","Use Amat (not Pmat) to extract diagonal fieldsplit blocks", "PCFieldSplitSetDiagUseAmat",jac->diag_use_amat,&jac->diag_use_amat,NULL);CHKERRQ(ierr);
12492686e3e9SMatthew G. Knepley   jac->offdiag_use_amat = pc->useAmat;
12502686e3e9SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_fieldsplit_off_diag_use_amat","Use Amat (not Pmat) to extract off-diagonal fieldsplit blocks", "PCFieldSplitSetOffDiagUseAmat",jac->offdiag_use_amat,&jac->offdiag_use_amat,NULL);CHKERRQ(ierr);
125100216929SDmitry Karpeev   /* FIXME: No programmatic equivalent to the following. */
1252c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1253c0adfefeSBarry Smith   if (stokes) {
1254c0adfefeSBarry Smith     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1255c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1256c0adfefeSBarry Smith   }
1257c0adfefeSBarry Smith 
12583b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
12593b224e63SBarry Smith   if (flg) {
12603b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
12613b224e63SBarry Smith   }
1262c30613efSMatthew Knepley   /* Only setup fields once */
1263c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1264d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1265d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
12666c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
12676c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1268d32f9abdSBarry Smith   }
1269c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1270c5929fdfSBarry Smith     ierr = PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1271c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
12720298fd71SBarry Smith     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,NULL);CHKERRQ(ierr);
127329f8a81cSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);CHKERRQ(ierr);
1274c096484dSStefano Zampini     ierr = PetscOptionsScalar("-pc_fieldsplit_schur_scale","Scale Schur complement","PCFieldSplitSetSchurScale",jac->schurscale,&jac->schurscale,NULL);CHKERRQ(ierr);
1275c5d2311dSJed Brown   }
12761b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
12770971522cSBarry Smith   PetscFunctionReturn(0);
12780971522cSBarry Smith }
12790971522cSBarry Smith 
12800971522cSBarry Smith /*------------------------------------------------------------------------------------*/
12810971522cSBarry Smith 
12821e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
12830971522cSBarry Smith {
128497bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
12850971522cSBarry Smith   PetscErrorCode    ierr;
12865a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
128769a612a9SBarry Smith   char              prefix[128];
12885d4c12cdSJungho Lee   PetscInt          i;
12890971522cSBarry Smith 
12900971522cSBarry Smith   PetscFunctionBegin;
12916c924f48SJed Brown   if (jac->splitdefined) {
12926c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
12936c924f48SJed Brown     PetscFunctionReturn(0);
12946c924f48SJed Brown   }
129551f519a2SBarry Smith   for (i=0; i<n; i++) {
1296e32f2f54SBarry 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);
1297e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
129851f519a2SBarry Smith   }
1299b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1300a04f6461SBarry Smith   if (splitname) {
1301db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1302a04f6461SBarry Smith   } else {
1303785e854fSJed Brown     ierr = PetscMalloc1(3,&ilink->splitname);CHKERRQ(ierr);
1304a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1305a04f6461SBarry Smith   }
13069df09d43SBarry Smith   ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + jac->nsplits : KSP_Solve_FS_0 + 4; /* Any split great than 4 gets logged in the 4th split */
1307785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields);CHKERRQ(ierr);
13085a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1309785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields_col);CHKERRQ(ierr);
13105d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
13112fa5cd67SKarl Rupp 
13125a9f2f41SSatish Balay   ilink->nfields = n;
13130298fd71SBarry Smith   ilink->next    = NULL;
1314ce94432eSBarry Smith   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1315422a814eSBarry Smith   ierr           = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
131620252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
13175a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
13189005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
131969a612a9SBarry Smith 
13208caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
13215a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
13220971522cSBarry Smith 
13230971522cSBarry Smith   if (!next) {
13245a9f2f41SSatish Balay     jac->head       = ilink;
13250298fd71SBarry Smith     ilink->previous = NULL;
13260971522cSBarry Smith   } else {
13270971522cSBarry Smith     while (next->next) {
13280971522cSBarry Smith       next = next->next;
13290971522cSBarry Smith     }
13305a9f2f41SSatish Balay     next->next      = ilink;
133151f519a2SBarry Smith     ilink->previous = next;
13320971522cSBarry Smith   }
13330971522cSBarry Smith   jac->nsplits++;
13340971522cSBarry Smith   PetscFunctionReturn(0);
13350971522cSBarry Smith }
13360971522cSBarry Smith 
13371e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1338e69d4d44SBarry Smith {
1339e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1340e69d4d44SBarry Smith   PetscErrorCode ierr;
1341e69d4d44SBarry Smith 
1342e69d4d44SBarry Smith   PetscFunctionBegin;
134334a3b412SBarry Smith   if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()");
1344785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
1345e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
13462fa5cd67SKarl Rupp 
1347e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
134813e0d083SBarry Smith   if (n) *n = jac->nsplits;
1349e69d4d44SBarry Smith   PetscFunctionReturn(0);
1350e69d4d44SBarry Smith }
13510971522cSBarry Smith 
13521e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
13530971522cSBarry Smith {
13540971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
13550971522cSBarry Smith   PetscErrorCode    ierr;
13560971522cSBarry Smith   PetscInt          cnt   = 0;
13575a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
13580971522cSBarry Smith 
13590971522cSBarry Smith   PetscFunctionBegin;
1360785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
13615a9f2f41SSatish Balay   while (ilink) {
13625a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
13635a9f2f41SSatish Balay     ilink            = ilink->next;
13640971522cSBarry Smith   }
13655d480477SMatthew 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);
136613e0d083SBarry Smith   if (n) *n = jac->nsplits;
13670971522cSBarry Smith   PetscFunctionReturn(0);
13680971522cSBarry Smith }
13690971522cSBarry Smith 
13706dbb499eSCian Wilson /*@C
13716dbb499eSCian Wilson     PCFieldSplitRestrictIS - Restricts the fieldsplit ISs to be within a given IS.
13726dbb499eSCian Wilson 
13736dbb499eSCian Wilson     Input Parameters:
13746dbb499eSCian Wilson +   pc  - the preconditioner context
13756dbb499eSCian Wilson +   is - the index set that defines the indices to which the fieldsplit is to be restricted
13766dbb499eSCian Wilson 
13776dbb499eSCian Wilson     Level: advanced
13786dbb499eSCian Wilson 
13796dbb499eSCian Wilson @*/
13806dbb499eSCian Wilson PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
13816dbb499eSCian Wilson {
13826dbb499eSCian Wilson   PetscErrorCode ierr;
13836dbb499eSCian Wilson 
13846dbb499eSCian Wilson   PetscFunctionBegin;
13856dbb499eSCian Wilson   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13866dbb499eSCian Wilson   PetscValidHeaderSpecific(isy,IS_CLASSID,2);
13870246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));CHKERRQ(ierr);
13886dbb499eSCian Wilson   PetscFunctionReturn(0);
13896dbb499eSCian Wilson }
13906dbb499eSCian Wilson 
13916dbb499eSCian Wilson 
13926dbb499eSCian Wilson static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
13936dbb499eSCian Wilson {
13946dbb499eSCian Wilson   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
13956dbb499eSCian Wilson   PetscErrorCode    ierr;
13966dbb499eSCian Wilson   PC_FieldSplitLink ilink = jac->head, next;
13976dbb499eSCian Wilson   PetscInt          localsize,size,sizez,i;
13986dbb499eSCian Wilson   const PetscInt    *ind, *indz;
13996dbb499eSCian Wilson   PetscInt          *indc, *indcz;
14006dbb499eSCian Wilson   PetscBool         flg;
14016dbb499eSCian Wilson 
14026dbb499eSCian Wilson   PetscFunctionBegin;
14036dbb499eSCian Wilson   ierr = ISGetLocalSize(isy,&localsize);CHKERRQ(ierr);
14046dbb499eSCian Wilson   ierr = MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));CHKERRQ(ierr);
14056dbb499eSCian Wilson   size -= localsize;
14066dbb499eSCian Wilson   while(ilink) {
14076dbb499eSCian Wilson     IS isrl,isr;
14081c7cfba8SBarry Smith     PC subpc;
14096dbb499eSCian Wilson     ierr          = ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);CHKERRQ(ierr);
14106dbb499eSCian Wilson     ierr          = ISGetLocalSize(isrl,&localsize);CHKERRQ(ierr);
14116dbb499eSCian Wilson     ierr          = PetscMalloc1(localsize,&indc);CHKERRQ(ierr);
14126dbb499eSCian Wilson     ierr          = ISGetIndices(isrl,&ind);CHKERRQ(ierr);
14136dbb499eSCian Wilson     ierr          = PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));CHKERRQ(ierr);
14146dbb499eSCian Wilson     ierr          = ISRestoreIndices(isrl,&ind);CHKERRQ(ierr);
14156dbb499eSCian Wilson     ierr          = ISDestroy(&isrl);CHKERRQ(ierr);
14166dbb499eSCian Wilson     for (i=0; i<localsize; i++) *(indc+i) += size;
14176dbb499eSCian Wilson     ierr          = ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);CHKERRQ(ierr);
14186dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
14196dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
14206dbb499eSCian Wilson     ilink->is     = isr;
14216dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
14226dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
14236dbb499eSCian Wilson     ilink->is_col = isr;
14246dbb499eSCian Wilson     ierr          = ISDestroy(&isr);CHKERRQ(ierr);
14256dbb499eSCian Wilson     ierr          = KSPGetPC(ilink->ksp, &subpc);CHKERRQ(ierr);
14266dbb499eSCian Wilson     ierr          = PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
14276dbb499eSCian Wilson     if(flg) {
14286dbb499eSCian Wilson       IS iszl,isz;
14296dbb499eSCian Wilson       MPI_Comm comm;
14306dbb499eSCian Wilson       ierr   = ISGetLocalSize(ilink->is,&localsize);CHKERRQ(ierr);
14316dbb499eSCian Wilson       comm   = PetscObjectComm((PetscObject)ilink->is);
14326dbb499eSCian Wilson       ierr   = ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);CHKERRQ(ierr);
14336dbb499eSCian Wilson       ierr   = MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
14346dbb499eSCian Wilson       sizez -= localsize;
14356dbb499eSCian Wilson       ierr   = ISGetLocalSize(iszl,&localsize);CHKERRQ(ierr);
14366dbb499eSCian Wilson       ierr   = PetscMalloc1(localsize,&indcz);CHKERRQ(ierr);
14376dbb499eSCian Wilson       ierr   = ISGetIndices(iszl,&indz);CHKERRQ(ierr);
14386dbb499eSCian Wilson       ierr   = PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));CHKERRQ(ierr);
14396dbb499eSCian Wilson       ierr   = ISRestoreIndices(iszl,&indz);CHKERRQ(ierr);
14406dbb499eSCian Wilson       ierr   = ISDestroy(&iszl);CHKERRQ(ierr);
14416dbb499eSCian Wilson       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
14426dbb499eSCian Wilson       ierr   = ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);CHKERRQ(ierr);
14436dbb499eSCian Wilson       ierr   = PCFieldSplitRestrictIS(subpc,isz);CHKERRQ(ierr);
14446dbb499eSCian Wilson       ierr   = ISDestroy(&isz);CHKERRQ(ierr);
14456dbb499eSCian Wilson     }
14466dbb499eSCian Wilson     next = ilink->next;
14476dbb499eSCian Wilson     ilink = next;
14486dbb499eSCian Wilson   }
14496dbb499eSCian Wilson   jac->isrestrict = PETSC_TRUE;
14506dbb499eSCian Wilson   PetscFunctionReturn(0);
14516dbb499eSCian Wilson }
14526dbb499eSCian Wilson 
14531e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1454704ba839SBarry Smith {
1455704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1456704ba839SBarry Smith   PetscErrorCode    ierr;
1457704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1458704ba839SBarry Smith   char              prefix[128];
1459704ba839SBarry Smith 
1460704ba839SBarry Smith   PetscFunctionBegin;
14616c924f48SJed Brown   if (jac->splitdefined) {
14626c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
14636c924f48SJed Brown     PetscFunctionReturn(0);
14646c924f48SJed Brown   }
1465b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1466a04f6461SBarry Smith   if (splitname) {
1467db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1468a04f6461SBarry Smith   } else {
1469785e854fSJed Brown     ierr = PetscMalloc1(8,&ilink->splitname);CHKERRQ(ierr);
1470b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1471a04f6461SBarry Smith   }
14729df09d43SBarry Smith   ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + jac->nsplits : KSP_Solve_FS_0 + 4; /* Any split great than 4 gets logged in the 4th split */
14731ab39975SBarry Smith   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1474b5787286SJed Brown   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1475b5787286SJed Brown   ilink->is     = is;
1476b5787286SJed Brown   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1477b5787286SJed Brown   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1478b5787286SJed Brown   ilink->is_col = is;
14790298fd71SBarry Smith   ilink->next   = NULL;
1480ce94432eSBarry Smith   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1481422a814eSBarry Smith   ierr          = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
148220252d06SBarry Smith   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1483704ba839SBarry Smith   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
14849005cf84SBarry Smith   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1485704ba839SBarry Smith 
14868caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1487704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1488704ba839SBarry Smith 
1489704ba839SBarry Smith   if (!next) {
1490704ba839SBarry Smith     jac->head       = ilink;
14910298fd71SBarry Smith     ilink->previous = NULL;
1492704ba839SBarry Smith   } else {
1493704ba839SBarry Smith     while (next->next) {
1494704ba839SBarry Smith       next = next->next;
1495704ba839SBarry Smith     }
1496704ba839SBarry Smith     next->next      = ilink;
1497704ba839SBarry Smith     ilink->previous = next;
1498704ba839SBarry Smith   }
1499704ba839SBarry Smith   jac->nsplits++;
1500704ba839SBarry Smith   PetscFunctionReturn(0);
1501704ba839SBarry Smith }
1502704ba839SBarry Smith 
150394ef8ddeSSatish Balay /*@C
15040971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
15050971522cSBarry Smith 
1506ad4df100SBarry Smith     Logically Collective on PC
15070971522cSBarry Smith 
15080971522cSBarry Smith     Input Parameters:
15090971522cSBarry Smith +   pc  - the preconditioner context
15100298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
15110971522cSBarry Smith .   n - the number of fields in this split
1512db4c96c1SJed Brown -   fields - the fields in this split
15130971522cSBarry Smith 
15140971522cSBarry Smith     Level: intermediate
15150971522cSBarry Smith 
1516*95452b02SPatrick Sanan     Notes:
1517*95452b02SPatrick Sanan     Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1518d32f9abdSBarry Smith 
15197287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1520d32f9abdSBarry 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
1521d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1522d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1523d32f9abdSBarry Smith 
1524db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1525db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1526db4c96c1SJed Brown 
15275d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
15285d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
15295d4c12cdSJungho Lee      available when this routine is called.
15305d4c12cdSJungho Lee 
1531d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
15320971522cSBarry Smith 
15330971522cSBarry Smith @*/
15345d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
15350971522cSBarry Smith {
15364ac538c5SBarry Smith   PetscErrorCode ierr;
15370971522cSBarry Smith 
15380971522cSBarry Smith   PetscFunctionBegin;
15390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1540db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1541ce94432eSBarry Smith   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1542db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
15430246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
15440971522cSBarry Smith   PetscFunctionReturn(0);
15450971522cSBarry Smith }
15460971522cSBarry Smith 
1547c84da90fSDmitry Karpeev /*@
1548c84da90fSDmitry Karpeev     PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1549c84da90fSDmitry Karpeev 
1550c84da90fSDmitry Karpeev     Logically Collective on PC
1551c84da90fSDmitry Karpeev 
1552c84da90fSDmitry Karpeev     Input Parameters:
1553c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1554c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1555c84da90fSDmitry Karpeev 
15569506b623SDmitry Karpeev     Options Database:
15579506b623SDmitry Karpeev .     -pc_fieldsplit_diag_use_amat
1558c84da90fSDmitry Karpeev 
1559c84da90fSDmitry Karpeev     Level: intermediate
1560c84da90fSDmitry Karpeev 
1561c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetDiagUseAmat(), PCFieldSplitSetOffDiagUseAmat(), PCFIELDSPLIT
1562c84da90fSDmitry Karpeev 
1563c84da90fSDmitry Karpeev @*/
1564c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1565c84da90fSDmitry Karpeev {
1566c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1567c84da90fSDmitry Karpeev   PetscBool      isfs;
1568c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1569c84da90fSDmitry Karpeev 
1570c84da90fSDmitry Karpeev   PetscFunctionBegin;
1571c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1572c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1573c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1574c84da90fSDmitry Karpeev   jac->diag_use_amat = flg;
1575c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1576c84da90fSDmitry Karpeev }
1577c84da90fSDmitry Karpeev 
1578c84da90fSDmitry Karpeev /*@
1579c84da90fSDmitry Karpeev     PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1580c84da90fSDmitry Karpeev 
1581c84da90fSDmitry Karpeev     Logically Collective on PC
1582c84da90fSDmitry Karpeev 
1583c84da90fSDmitry Karpeev     Input Parameters:
1584c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1585c84da90fSDmitry Karpeev 
1586c84da90fSDmitry Karpeev     Output Parameters:
1587c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1588c84da90fSDmitry Karpeev 
1589c84da90fSDmitry Karpeev 
1590c84da90fSDmitry Karpeev     Level: intermediate
1591c84da90fSDmitry Karpeev 
1592c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetDiagUseAmat(), PCFieldSplitGetOffDiagUseAmat(), PCFIELDSPLIT
1593c84da90fSDmitry Karpeev 
1594c84da90fSDmitry Karpeev @*/
1595c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1596c84da90fSDmitry Karpeev {
1597c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1598c84da90fSDmitry Karpeev   PetscBool      isfs;
1599c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1600c84da90fSDmitry Karpeev 
1601c84da90fSDmitry Karpeev   PetscFunctionBegin;
1602c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1603c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1604c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1605c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1606c84da90fSDmitry Karpeev   *flg = jac->diag_use_amat;
1607c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1608c84da90fSDmitry Karpeev }
1609c84da90fSDmitry Karpeev 
1610c84da90fSDmitry Karpeev /*@
1611c84da90fSDmitry Karpeev     PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1612c84da90fSDmitry Karpeev 
1613c84da90fSDmitry Karpeev     Logically Collective on PC
1614c84da90fSDmitry Karpeev 
1615c84da90fSDmitry Karpeev     Input Parameters:
1616c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1617c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1618c84da90fSDmitry Karpeev 
16199506b623SDmitry Karpeev     Options Database:
16209506b623SDmitry Karpeev .     -pc_fieldsplit_off_diag_use_amat
1621c84da90fSDmitry Karpeev 
1622c84da90fSDmitry Karpeev     Level: intermediate
1623c84da90fSDmitry Karpeev 
1624c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetOffDiagUseAmat(), PCFieldSplitSetDiagUseAmat(), PCFIELDSPLIT
1625c84da90fSDmitry Karpeev 
1626c84da90fSDmitry Karpeev @*/
1627c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1628c84da90fSDmitry Karpeev {
1629c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1630c84da90fSDmitry Karpeev   PetscBool      isfs;
1631c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1632c84da90fSDmitry Karpeev 
1633c84da90fSDmitry Karpeev   PetscFunctionBegin;
1634c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1635c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1636c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1637c84da90fSDmitry Karpeev   jac->offdiag_use_amat = flg;
1638c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1639c84da90fSDmitry Karpeev }
1640c84da90fSDmitry Karpeev 
1641c84da90fSDmitry Karpeev /*@
1642c84da90fSDmitry Karpeev     PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1643c84da90fSDmitry Karpeev 
1644c84da90fSDmitry Karpeev     Logically Collective on PC
1645c84da90fSDmitry Karpeev 
1646c84da90fSDmitry Karpeev     Input Parameters:
1647c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1648c84da90fSDmitry Karpeev 
1649c84da90fSDmitry Karpeev     Output Parameters:
1650c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1651c84da90fSDmitry Karpeev 
1652c84da90fSDmitry Karpeev 
1653c84da90fSDmitry Karpeev     Level: intermediate
1654c84da90fSDmitry Karpeev 
1655c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetOffDiagUseAmat(), PCFieldSplitGetDiagUseAmat(), PCFIELDSPLIT
1656c84da90fSDmitry Karpeev 
1657c84da90fSDmitry Karpeev @*/
1658c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1659c84da90fSDmitry Karpeev {
1660c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1661c84da90fSDmitry Karpeev   PetscBool      isfs;
1662c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1663c84da90fSDmitry Karpeev 
1664c84da90fSDmitry Karpeev   PetscFunctionBegin;
1665c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1666c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1667c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1668c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1669c84da90fSDmitry Karpeev   *flg = jac->offdiag_use_amat;
1670c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1671c84da90fSDmitry Karpeev }
1672c84da90fSDmitry Karpeev 
1673c84da90fSDmitry Karpeev 
1674c84da90fSDmitry Karpeev 
1675218d4943SBarry Smith /*@C
1676704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1677704ba839SBarry Smith 
1678ad4df100SBarry Smith     Logically Collective on PC
1679704ba839SBarry Smith 
1680704ba839SBarry Smith     Input Parameters:
1681704ba839SBarry Smith +   pc  - the preconditioner context
16820298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
1683db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1684704ba839SBarry Smith 
1685d32f9abdSBarry Smith 
1686a6ffb8dbSJed Brown     Notes:
1687a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1688a6ffb8dbSJed Brown 
1689db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1690db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1691d32f9abdSBarry Smith 
1692704ba839SBarry Smith     Level: intermediate
1693704ba839SBarry Smith 
1694704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1695704ba839SBarry Smith 
1696704ba839SBarry Smith @*/
16977087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1698704ba839SBarry Smith {
16994ac538c5SBarry Smith   PetscErrorCode ierr;
1700704ba839SBarry Smith 
1701704ba839SBarry Smith   PetscFunctionBegin;
17020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17037b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1704db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
1705ccc738f7SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1706704ba839SBarry Smith   PetscFunctionReturn(0);
1707704ba839SBarry Smith }
1708704ba839SBarry Smith 
170994ef8ddeSSatish Balay /*@C
171057a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
171157a9adfeSMatthew G Knepley 
171257a9adfeSMatthew G Knepley     Logically Collective on PC
171357a9adfeSMatthew G Knepley 
171457a9adfeSMatthew G Knepley     Input Parameters:
171557a9adfeSMatthew G Knepley +   pc  - the preconditioner context
171657a9adfeSMatthew G Knepley -   splitname - name of this split
171757a9adfeSMatthew G Knepley 
171857a9adfeSMatthew G Knepley     Output Parameter:
17190298fd71SBarry Smith -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
172057a9adfeSMatthew G Knepley 
172157a9adfeSMatthew G Knepley     Level: intermediate
172257a9adfeSMatthew G Knepley 
172357a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
172457a9adfeSMatthew G Knepley 
172557a9adfeSMatthew G Knepley @*/
172657a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
172757a9adfeSMatthew G Knepley {
172857a9adfeSMatthew G Knepley   PetscErrorCode ierr;
172957a9adfeSMatthew G Knepley 
173057a9adfeSMatthew G Knepley   PetscFunctionBegin;
173157a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
173257a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
173357a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
173457a9adfeSMatthew G Knepley   {
173557a9adfeSMatthew G Knepley     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
173657a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
173757a9adfeSMatthew G Knepley     PetscBool         found;
173857a9adfeSMatthew G Knepley 
17390298fd71SBarry Smith     *is = NULL;
174057a9adfeSMatthew G Knepley     while (ilink) {
174157a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
174257a9adfeSMatthew G Knepley       if (found) {
174357a9adfeSMatthew G Knepley         *is = ilink->is;
174457a9adfeSMatthew G Knepley         break;
174557a9adfeSMatthew G Knepley       }
174657a9adfeSMatthew G Knepley       ilink = ilink->next;
174757a9adfeSMatthew G Knepley     }
174857a9adfeSMatthew G Knepley   }
174957a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
175057a9adfeSMatthew G Knepley }
175157a9adfeSMatthew G Knepley 
175251f519a2SBarry Smith /*@
175351f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
175451f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
175551f519a2SBarry Smith 
1756ad4df100SBarry Smith     Logically Collective on PC
175751f519a2SBarry Smith 
175851f519a2SBarry Smith     Input Parameters:
175951f519a2SBarry Smith +   pc  - the preconditioner context
176051f519a2SBarry Smith -   bs - the block size
176151f519a2SBarry Smith 
176251f519a2SBarry Smith     Level: intermediate
176351f519a2SBarry Smith 
176451f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
176551f519a2SBarry Smith 
176651f519a2SBarry Smith @*/
17677087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
176851f519a2SBarry Smith {
17694ac538c5SBarry Smith   PetscErrorCode ierr;
177051f519a2SBarry Smith 
177151f519a2SBarry Smith   PetscFunctionBegin;
17720700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1773c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
17744ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
177551f519a2SBarry Smith   PetscFunctionReturn(0);
177651f519a2SBarry Smith }
177751f519a2SBarry Smith 
17780971522cSBarry Smith /*@C
177969a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
17800971522cSBarry Smith 
178169a612a9SBarry Smith    Collective on KSP
17820971522cSBarry Smith 
17830971522cSBarry Smith    Input Parameter:
17840971522cSBarry Smith .  pc - the preconditioner context
17850971522cSBarry Smith 
17860971522cSBarry Smith    Output Parameters:
178713e0d083SBarry Smith +  n - the number of splits
17883111de3cSBarry Smith -  subksp - the array of KSP contexts
17890971522cSBarry Smith 
17900971522cSBarry Smith    Note:
17919a4f7e47SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
1792d32f9abdSBarry Smith    (not the KSP just the array that contains them).
17930971522cSBarry Smith 
179469a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
17950971522cSBarry Smith 
1796196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
17972bf68e3eSBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
1798196cc216SBarry Smith       KSP array must be.
1799196cc216SBarry Smith 
1800196cc216SBarry Smith 
18010971522cSBarry Smith    Level: advanced
18020971522cSBarry Smith 
18030971522cSBarry Smith .seealso: PCFIELDSPLIT
18040971522cSBarry Smith @*/
18057087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
18060971522cSBarry Smith {
18074ac538c5SBarry Smith   PetscErrorCode ierr;
18080971522cSBarry Smith 
18090971522cSBarry Smith   PetscFunctionBegin;
18100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181113e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
18124ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
18130971522cSBarry Smith   PetscFunctionReturn(0);
18140971522cSBarry Smith }
18150971522cSBarry Smith 
1816e69d4d44SBarry Smith /*@
181753f2277eSBarry Smith     PCFieldSplitSetSchurPre -  Indicates what operator is used to construct the preconditioner for the Schur complement.
1818a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1819e69d4d44SBarry Smith 
1820e69d4d44SBarry Smith     Collective on PC
1821e69d4d44SBarry Smith 
1822e69d4d44SBarry Smith     Input Parameters:
1823e69d4d44SBarry Smith +   pc      - the preconditioner context
18246fdd48a9SBarry Smith .   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_SCHUR_PRE_USER
18256fdd48a9SBarry Smith               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
18260298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
1827084e4875SJed Brown 
1828e69d4d44SBarry Smith     Options Database:
18296fdd48a9SBarry Smith .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments
1830e69d4d44SBarry Smith 
1831fd1303e9SJungho Lee     Notes:
1832fd1303e9SJungho Lee $    If ptype is
1833a6a584a2SBarry Smith $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
183453f2277eSBarry Smith $             matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
1835a6a584a2SBarry Smith $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
1836a6a584a2SBarry Smith $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
1837fd1303e9SJungho Lee $             preconditioner
18386fdd48a9SBarry Smith $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
18396fdd48a9SBarry Smith $             to this function).
1840a6a584a2SBarry Smith $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
18411d71051fSDmitry Karpeev $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
18426fdd48a9SBarry Smith $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
18436fdd48a9SBarry Smith $        full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PFIELDSPLIT (this is expensive)
18446fdd48a9SBarry Smith $             useful mostly as a test that the Schur complement approach can work for your problem
1845fd1303e9SJungho Lee 
1846e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1847fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1848a7476a74SDmitry Karpeev     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
1849fd1303e9SJungho Lee 
1850e69d4d44SBarry Smith     Level: intermediate
1851e69d4d44SBarry Smith 
18521d71051fSDmitry Karpeev .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
18531d71051fSDmitry Karpeev           MatSchurComplementSetAinvType(), PCLSC
1854e69d4d44SBarry Smith 
1855e69d4d44SBarry Smith @*/
185629f8a81cSJed Brown PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1857e69d4d44SBarry Smith {
18584ac538c5SBarry Smith   PetscErrorCode ierr;
1859e69d4d44SBarry Smith 
1860e69d4d44SBarry Smith   PetscFunctionBegin;
18610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186229f8a81cSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1863e69d4d44SBarry Smith   PetscFunctionReturn(0);
1864e69d4d44SBarry Smith }
186529f8a81cSJed Brown PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */
1866e69d4d44SBarry Smith 
186737a82bf0SJed Brown /*@
186837a82bf0SJed Brown     PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
186929f8a81cSJed Brown     preconditioned.  See PCFieldSplitSetSchurPre() for details.
187037a82bf0SJed Brown 
187137a82bf0SJed Brown     Logically Collective on PC
187237a82bf0SJed Brown 
187337a82bf0SJed Brown     Input Parameters:
187437a82bf0SJed Brown .   pc      - the preconditioner context
187537a82bf0SJed Brown 
187637a82bf0SJed Brown     Output Parameters:
187737a82bf0SJed Brown +   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11, PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
187837a82bf0SJed Brown -   userpre - matrix to use for preconditioning (with PC_FIELDSPLIT_PRE_USER), or NULL
187937a82bf0SJed Brown 
188037a82bf0SJed Brown     Level: intermediate
188137a82bf0SJed Brown 
188229f8a81cSJed Brown .seealso: PCFieldSplitSetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
188337a82bf0SJed Brown 
188437a82bf0SJed Brown @*/
188537a82bf0SJed Brown PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
188637a82bf0SJed Brown {
188737a82bf0SJed Brown   PetscErrorCode ierr;
188837a82bf0SJed Brown 
188937a82bf0SJed Brown   PetscFunctionBegin;
189037a82bf0SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
189137a82bf0SJed Brown   ierr = PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));CHKERRQ(ierr);
1892e69d4d44SBarry Smith   PetscFunctionReturn(0);
1893e69d4d44SBarry Smith }
1894e69d4d44SBarry Smith 
189545e7fc46SDmitry Karpeev /*@
1896470b340bSDmitry Karpeev     PCFieldSplitSchurGetS -  extract the MatSchurComplement object used by this PC in case it needs to be configured separately
189745e7fc46SDmitry Karpeev 
189845e7fc46SDmitry Karpeev     Not collective
189945e7fc46SDmitry Karpeev 
190045e7fc46SDmitry Karpeev     Input Parameter:
190145e7fc46SDmitry Karpeev .   pc      - the preconditioner context
190245e7fc46SDmitry Karpeev 
190345e7fc46SDmitry Karpeev     Output Parameter:
1904470b340bSDmitry Karpeev .   S       - the Schur complement matrix
190545e7fc46SDmitry Karpeev 
1906470b340bSDmitry Karpeev     Notes:
1907470b340bSDmitry Karpeev     This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().
190845e7fc46SDmitry Karpeev 
190945e7fc46SDmitry Karpeev     Level: advanced
191045e7fc46SDmitry Karpeev 
191129f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurRestoreS()
191245e7fc46SDmitry Karpeev 
191345e7fc46SDmitry Karpeev @*/
1914470b340bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
191545e7fc46SDmitry Karpeev {
191645e7fc46SDmitry Karpeev   PetscErrorCode ierr;
191745e7fc46SDmitry Karpeev   const char*    t;
191845e7fc46SDmitry Karpeev   PetscBool      isfs;
191945e7fc46SDmitry Karpeev   PC_FieldSplit  *jac;
192045e7fc46SDmitry Karpeev 
192145e7fc46SDmitry Karpeev   PetscFunctionBegin;
192245e7fc46SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
192345e7fc46SDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
192445e7fc46SDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
192545e7fc46SDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
192645e7fc46SDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
1927470b340bSDmitry Karpeev   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1928470b340bSDmitry Karpeev   if (S) *S = jac->schur;
192945e7fc46SDmitry Karpeev   PetscFunctionReturn(0);
193045e7fc46SDmitry Karpeev }
193145e7fc46SDmitry Karpeev 
1932470b340bSDmitry Karpeev /*@
1933470b340bSDmitry Karpeev     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC
1934470b340bSDmitry Karpeev 
1935470b340bSDmitry Karpeev     Not collective
1936470b340bSDmitry Karpeev 
1937470b340bSDmitry Karpeev     Input Parameters:
1938470b340bSDmitry Karpeev +   pc      - the preconditioner context
1939470b340bSDmitry Karpeev .   S       - the Schur complement matrix
1940470b340bSDmitry Karpeev 
1941470b340bSDmitry Karpeev     Level: advanced
1942470b340bSDmitry Karpeev 
194329f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurGetS()
1944470b340bSDmitry Karpeev 
1945470b340bSDmitry Karpeev @*/
1946bca69d2bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
1947470b340bSDmitry Karpeev {
1948470b340bSDmitry Karpeev   PetscErrorCode ierr;
1949470b340bSDmitry Karpeev   const char*    t;
1950470b340bSDmitry Karpeev   PetscBool      isfs;
1951470b340bSDmitry Karpeev   PC_FieldSplit  *jac;
1952470b340bSDmitry Karpeev 
1953470b340bSDmitry Karpeev   PetscFunctionBegin;
1954470b340bSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1955470b340bSDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
1956470b340bSDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1957470b340bSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1958470b340bSDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
1959470b340bSDmitry Karpeev   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1960bca69d2bSDmitry Karpeev   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
1961470b340bSDmitry Karpeev   PetscFunctionReturn(0);
1962470b340bSDmitry Karpeev }
1963470b340bSDmitry Karpeev 
1964470b340bSDmitry Karpeev 
196529f8a81cSJed Brown static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1966e69d4d44SBarry Smith {
1967e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1968084e4875SJed Brown   PetscErrorCode ierr;
1969e69d4d44SBarry Smith 
1970e69d4d44SBarry Smith   PetscFunctionBegin;
1971084e4875SJed Brown   jac->schurpre = ptype;
1972a7476a74SDmitry Karpeev   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
19736bf464f9SBarry Smith     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1974084e4875SJed Brown     jac->schur_user = pre;
1975084e4875SJed Brown     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1976084e4875SJed Brown   }
1977e69d4d44SBarry Smith   PetscFunctionReturn(0);
1978e69d4d44SBarry Smith }
1979e69d4d44SBarry Smith 
198037a82bf0SJed Brown static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
198137a82bf0SJed Brown {
198237a82bf0SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
198337a82bf0SJed Brown 
198437a82bf0SJed Brown   PetscFunctionBegin;
198537a82bf0SJed Brown   *ptype = jac->schurpre;
198637a82bf0SJed Brown   *pre   = jac->schur_user;
198737a82bf0SJed Brown   PetscFunctionReturn(0);
198837a82bf0SJed Brown }
198937a82bf0SJed Brown 
1990ab1df9f5SJed Brown /*@
19910ffb0e17SBarry Smith     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain in the preconditioner
1992ab1df9f5SJed Brown 
1993ab1df9f5SJed Brown     Collective on PC
1994ab1df9f5SJed Brown 
1995ab1df9f5SJed Brown     Input Parameters:
1996ab1df9f5SJed Brown +   pc  - the preconditioner context
1997c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1998ab1df9f5SJed Brown 
1999ab1df9f5SJed Brown     Options Database:
2000c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
2001ab1df9f5SJed Brown 
2002ab1df9f5SJed Brown 
2003ab1df9f5SJed Brown     Level: intermediate
2004ab1df9f5SJed Brown 
2005ab1df9f5SJed Brown     Notes:
2006ab1df9f5SJed Brown     The FULL factorization is
2007ab1df9f5SJed Brown 
20080ffb0e17SBarry Smith $   (A   B)  = (1       0) (A   0) (1  Ainv*B)  = L D U
20090ffb0e17SBarry Smith $   (C   E)    (C*Ainv  1) (0   S) (0     1  )
2010ab1df9f5SJed Brown 
20110ffb0e17SBarry Smith     where S = E - 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,
2012c096484dSStefano Zampini     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). Sign flipping of S can be turned off with PCFieldSplitSetSchurScale().
2013ab1df9f5SJed Brown 
20140ffb0e17SBarry Smith $    If A and S are solved exactly
20150ffb0e17SBarry Smith $      *) FULL factorization is a direct solver.
20160ffb0e17SBarry Smith $      *) The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial of degree 2, so KSPGMRES converges in 2 iterations.
20170ffb0e17SBarry Smith $      *) With DIAG, the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.
2018ab1df9f5SJed Brown 
20190ffb0e17SBarry Smith     If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
20200ffb0e17SBarry 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.
20210ffb0e17SBarry Smith 
20220ffb0e17SBarry Smith     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES.
20230ffb0e17SBarry Smith 
20240ffb0e17SBarry Smith     Note that a flexible method like KSPFGMRES or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used to solve with A or S).
2025ab1df9f5SJed Brown 
2026ab1df9f5SJed Brown     References:
202796a0c994SBarry Smith +   1. - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
202896a0c994SBarry Smith -   2. - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).
2029ab1df9f5SJed Brown 
2030c096484dSStefano Zampini .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCFieldSplitSetSchurScale()
2031ab1df9f5SJed Brown @*/
2032c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2033ab1df9f5SJed Brown {
2034ab1df9f5SJed Brown   PetscErrorCode ierr;
2035ab1df9f5SJed Brown 
2036ab1df9f5SJed Brown   PetscFunctionBegin;
2037ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2038c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
2039ab1df9f5SJed Brown   PetscFunctionReturn(0);
2040ab1df9f5SJed Brown }
2041ab1df9f5SJed Brown 
20421e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2043ab1df9f5SJed Brown {
2044ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2045ab1df9f5SJed Brown 
2046ab1df9f5SJed Brown   PetscFunctionBegin;
2047ab1df9f5SJed Brown   jac->schurfactorization = ftype;
2048ab1df9f5SJed Brown   PetscFunctionReturn(0);
2049ab1df9f5SJed Brown }
2050ab1df9f5SJed Brown 
2051c096484dSStefano Zampini /*@
2052c096484dSStefano Zampini     PCFieldSplitSetSchurScale -  Controls the sign flip of S for PC_FIELDSPLIT_SCHUR_FACT_DIAG.
2053c096484dSStefano Zampini 
2054c096484dSStefano Zampini     Collective on PC
2055c096484dSStefano Zampini 
2056c096484dSStefano Zampini     Input Parameters:
2057c096484dSStefano Zampini +   pc    - the preconditioner context
2058c096484dSStefano Zampini -   scale - scaling factor for the Schur complement
2059c096484dSStefano Zampini 
2060c096484dSStefano Zampini     Options Database:
2061c096484dSStefano Zampini .     -pc_fieldsplit_schur_scale - default is -1.0
2062c096484dSStefano Zampini 
2063c096484dSStefano Zampini     Level: intermediate
2064c096484dSStefano Zampini 
2065c096484dSStefano Zampini .seealso: PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurFactType, PCFieldSplitSetSchurScale()
2066c096484dSStefano Zampini @*/
2067c096484dSStefano Zampini PetscErrorCode PCFieldSplitSetSchurScale(PC pc,PetscScalar scale)
2068c096484dSStefano Zampini {
2069c096484dSStefano Zampini   PetscErrorCode ierr;
2070c096484dSStefano Zampini 
2071c096484dSStefano Zampini   PetscFunctionBegin;
2072c096484dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2073c096484dSStefano Zampini   PetscValidLogicalCollectiveScalar(pc,scale,2);
2074c096484dSStefano Zampini   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurScale_C",(PC,PetscScalar),(pc,scale));CHKERRQ(ierr);
2075c096484dSStefano Zampini   PetscFunctionReturn(0);
2076c096484dSStefano Zampini }
2077c096484dSStefano Zampini 
2078c096484dSStefano Zampini static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc,PetscScalar scale)
2079c096484dSStefano Zampini {
2080c096484dSStefano Zampini   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2081c096484dSStefano Zampini 
2082c096484dSStefano Zampini   PetscFunctionBegin;
2083c096484dSStefano Zampini   jac->schurscale = scale;
2084c096484dSStefano Zampini   PetscFunctionReturn(0);
2085c096484dSStefano Zampini }
2086c096484dSStefano Zampini 
208730ad9308SMatthew Knepley /*@C
20888c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
208930ad9308SMatthew Knepley 
209030ad9308SMatthew Knepley    Collective on KSP
209130ad9308SMatthew Knepley 
209230ad9308SMatthew Knepley    Input Parameter:
209330ad9308SMatthew Knepley .  pc - the preconditioner context
209430ad9308SMatthew Knepley 
209530ad9308SMatthew Knepley    Output Parameters:
2096a04f6461SBarry Smith +  A00 - the (0,0) block
2097a04f6461SBarry Smith .  A01 - the (0,1) block
2098a04f6461SBarry Smith .  A10 - the (1,0) block
2099a04f6461SBarry Smith -  A11 - the (1,1) block
210030ad9308SMatthew Knepley 
210130ad9308SMatthew Knepley    Level: advanced
210230ad9308SMatthew Knepley 
210330ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
210430ad9308SMatthew Knepley @*/
2105a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
210630ad9308SMatthew Knepley {
210730ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
210830ad9308SMatthew Knepley 
210930ad9308SMatthew Knepley   PetscFunctionBegin;
21100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2111ce94432eSBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2112a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
2113a04f6461SBarry Smith   if (A01) *A01 = jac->B;
2114a04f6461SBarry Smith   if (A10) *A10 = jac->C;
2115a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
211630ad9308SMatthew Knepley   PetscFunctionReturn(0);
211730ad9308SMatthew Knepley }
211830ad9308SMatthew Knepley 
21191e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
212079416396SBarry Smith {
212179416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2122e69d4d44SBarry Smith   PetscErrorCode ierr;
212379416396SBarry Smith 
212479416396SBarry Smith   PetscFunctionBegin;
212579416396SBarry Smith   jac->type = type;
21263b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
21273b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
21283b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
21292fa5cd67SKarl Rupp 
2130bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
213129f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);CHKERRQ(ierr);
213237a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);CHKERRQ(ierr);
2133bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
2134c096484dSStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",PCFieldSplitSetSchurScale_FieldSplit);CHKERRQ(ierr);
2135e69d4d44SBarry Smith 
21363b224e63SBarry Smith   } else {
21373b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
21383b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
21392fa5cd67SKarl Rupp 
2140bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
214129f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);CHKERRQ(ierr);
214237a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);CHKERRQ(ierr);
2143bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
2144c096484dSStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",0);CHKERRQ(ierr);
21453b224e63SBarry Smith   }
214679416396SBarry Smith   PetscFunctionReturn(0);
214779416396SBarry Smith }
214879416396SBarry Smith 
21491e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
215051f519a2SBarry Smith {
215151f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
215251f519a2SBarry Smith 
215351f519a2SBarry Smith   PetscFunctionBegin;
2154ce94432eSBarry Smith   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2155ce94432eSBarry Smith   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
215651f519a2SBarry Smith   jac->bs = bs;
215751f519a2SBarry Smith   PetscFunctionReturn(0);
215851f519a2SBarry Smith }
215951f519a2SBarry Smith 
2160bc08b0f1SBarry Smith /*@
216179416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
216279416396SBarry Smith 
216379416396SBarry Smith    Collective on PC
216479416396SBarry Smith 
216579416396SBarry Smith    Input Parameter:
216679416396SBarry Smith .  pc - the preconditioner context
216781540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
216879416396SBarry Smith 
216979416396SBarry Smith    Options Database Key:
2170a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
217179416396SBarry Smith 
2172b02e2d75SMatthew G Knepley    Level: Intermediate
217379416396SBarry Smith 
217479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
217579416396SBarry Smith 
217679416396SBarry Smith .seealso: PCCompositeSetType()
217779416396SBarry Smith 
217879416396SBarry Smith @*/
21797087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
218079416396SBarry Smith {
21814ac538c5SBarry Smith   PetscErrorCode ierr;
218279416396SBarry Smith 
218379416396SBarry Smith   PetscFunctionBegin;
21840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
21854ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
218679416396SBarry Smith   PetscFunctionReturn(0);
218779416396SBarry Smith }
218879416396SBarry Smith 
2189b02e2d75SMatthew G Knepley /*@
2190b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
2191b02e2d75SMatthew G Knepley 
2192b02e2d75SMatthew G Knepley   Not collective
2193b02e2d75SMatthew G Knepley 
2194b02e2d75SMatthew G Knepley   Input Parameter:
2195b02e2d75SMatthew G Knepley . pc - the preconditioner context
2196b02e2d75SMatthew G Knepley 
2197b02e2d75SMatthew G Knepley   Output Parameter:
2198b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
2199b02e2d75SMatthew G Knepley 
2200b02e2d75SMatthew G Knepley   Level: Intermediate
2201b02e2d75SMatthew G Knepley 
2202b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2203b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
2204b02e2d75SMatthew G Knepley @*/
2205b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2206b02e2d75SMatthew G Knepley {
2207b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2208b02e2d75SMatthew G Knepley 
2209b02e2d75SMatthew G Knepley   PetscFunctionBegin;
2210b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2211b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
2212b02e2d75SMatthew G Knepley   *type = jac->type;
2213b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
2214b02e2d75SMatthew G Knepley }
2215b02e2d75SMatthew G Knepley 
22164ab8060aSDmitry Karpeev /*@
22174ab8060aSDmitry Karpeev    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
22184ab8060aSDmitry Karpeev 
22194ab8060aSDmitry Karpeev    Logically Collective
22204ab8060aSDmitry Karpeev 
22214ab8060aSDmitry Karpeev    Input Parameters:
22224ab8060aSDmitry Karpeev +  pc   - the preconditioner context
22234ab8060aSDmitry Karpeev -  flg  - boolean indicating whether to use field splits defined by the DM
22244ab8060aSDmitry Karpeev 
22254ab8060aSDmitry Karpeev    Options Database Key:
22264ab8060aSDmitry Karpeev .  -pc_fieldsplit_dm_splits
22274ab8060aSDmitry Karpeev 
22284ab8060aSDmitry Karpeev    Level: Intermediate
22294ab8060aSDmitry Karpeev 
22304ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
22314ab8060aSDmitry Karpeev 
22324ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits()
22334ab8060aSDmitry Karpeev 
22344ab8060aSDmitry Karpeev @*/
22354ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
22364ab8060aSDmitry Karpeev {
22374ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
22384ab8060aSDmitry Karpeev   PetscBool      isfs;
22394ab8060aSDmitry Karpeev   PetscErrorCode ierr;
22404ab8060aSDmitry Karpeev 
22414ab8060aSDmitry Karpeev   PetscFunctionBegin;
22424ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22434ab8060aSDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
22444ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
22454ab8060aSDmitry Karpeev   if (isfs) {
22464ab8060aSDmitry Karpeev     jac->dm_splits = flg;
22474ab8060aSDmitry Karpeev   }
22484ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
22494ab8060aSDmitry Karpeev }
22504ab8060aSDmitry Karpeev 
22514ab8060aSDmitry Karpeev 
22524ab8060aSDmitry Karpeev /*@
22534ab8060aSDmitry Karpeev    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
22544ab8060aSDmitry Karpeev 
22554ab8060aSDmitry Karpeev    Logically Collective
22564ab8060aSDmitry Karpeev 
22574ab8060aSDmitry Karpeev    Input Parameter:
22584ab8060aSDmitry Karpeev .  pc   - the preconditioner context
22594ab8060aSDmitry Karpeev 
22604ab8060aSDmitry Karpeev    Output Parameter:
22614ab8060aSDmitry Karpeev .  flg  - boolean indicating whether to use field splits defined by the DM
22624ab8060aSDmitry Karpeev 
22634ab8060aSDmitry Karpeev    Level: Intermediate
22644ab8060aSDmitry Karpeev 
22654ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
22664ab8060aSDmitry Karpeev 
22674ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits()
22684ab8060aSDmitry Karpeev 
22694ab8060aSDmitry Karpeev @*/
22704ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
22714ab8060aSDmitry Karpeev {
22724ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
22734ab8060aSDmitry Karpeev   PetscBool      isfs;
22744ab8060aSDmitry Karpeev   PetscErrorCode ierr;
22754ab8060aSDmitry Karpeev 
22764ab8060aSDmitry Karpeev   PetscFunctionBegin;
22774ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22784ab8060aSDmitry Karpeev   PetscValidPointer(flg,2);
22794ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
22804ab8060aSDmitry Karpeev   if (isfs) {
22814ab8060aSDmitry Karpeev     if(flg) *flg = jac->dm_splits;
22824ab8060aSDmitry Karpeev   }
22834ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
22844ab8060aSDmitry Karpeev }
22854ab8060aSDmitry Karpeev 
22860971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
22870971522cSBarry Smith /*MC
2288a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2289a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
22900971522cSBarry Smith 
2291edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
2292edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
22930971522cSBarry Smith 
2294a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
229569a612a9SBarry Smith          and set the options directly on the resulting KSP object
22960971522cSBarry Smith 
22970971522cSBarry Smith    Level: intermediate
22980971522cSBarry Smith 
229979416396SBarry Smith    Options Database Keys:
230081540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
230181540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
230281540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
230381540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
23040f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
2305a6a584a2SBarry Smith .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see PCFieldSplitSetSchurPre()
2306435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
230779416396SBarry Smith 
23085d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
23095d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
23105d4c12cdSJungho Lee 
2311c8a0d604SMatthew G Knepley    Notes:
2312c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
2313d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
2314d32f9abdSBarry Smith 
2315d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
2316d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
2317d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
2318d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
2319d32f9abdSBarry Smith 
2320c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
2321c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
2322c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
232313b05affSJed Brown $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2324c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
232513b05affSJed Brown      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
232613b05affSJed Brown $              S = A11 - A10 ksp(A00) A01
232713b05affSJed Brown      which is usually dense and not stored explicitly.  The action of ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
2328c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2329c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2330a6a584a2SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.
23311d71051fSDmitry Karpeev 
2332a7476a74SDmitry Karpeev      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
23335668aaf4SBarry Smith      diag gives
2334c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
2335c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
2336c096484dSStefano Zampini      note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. For SPD matrices J, the sign flip
2337c096484dSStefano Zampini      can be turned off with PCFieldSplitSetSchurScale() or by command line -pc_fieldsplit_schur_scale 1.0. The lower factorization is the inverse of
2338c8a0d604SMatthew G Knepley $              (  A00   0 )
2339c8a0d604SMatthew G Knepley $              (  A10   S )
2340c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2341c8a0d604SMatthew G Knepley $              ( A00 A01 )
2342c8a0d604SMatthew G Knepley $              (  0   S  )
2343c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
2344e69d4d44SBarry Smith 
2345edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
2346edf189efSBarry Smith      is used automatically for a second block.
2347edf189efSBarry Smith 
2348ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
2349ff218e97SBarry Smith      Generally it should be used with the AIJ format.
2350ff218e97SBarry Smith 
2351ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
2352ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
2353ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
23540716a85fSBarry Smith 
2355a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
23560971522cSBarry Smith 
23573bc631e0SBarry Smith    There is a nice discussion of block preconditioners in
23583bc631e0SBarry Smith 
23593bc631e0SBarry Smith [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
23603bc631e0SBarry Smith        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
23613bc631e0SBarry Smith        http://chess.cs.umd.edu/~elman/papers/tax.pdf
23623bc631e0SBarry Smith 
23634d308398SBarry Smith    The Constrained Pressure Preconditioner (CPR) can be implemented using PCCOMPOSITE with PCGALERKIN. CPR first solves an R A P subsystem, updates the
23644d308398SBarry Smith    residual on all variables (PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE)), and then applies a simple ILU like preconditioner on all the variables.
2365a6a584a2SBarry Smith 
23667e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
23671d71051fSDmitry Karpeev            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
2368c096484dSStefano Zampini           MatSchurComplementSetAinvType(), PCFieldSplitSetSchurScale()
23690971522cSBarry Smith M*/
23700971522cSBarry Smith 
23718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
23720971522cSBarry Smith {
23730971522cSBarry Smith   PetscErrorCode ierr;
23740971522cSBarry Smith   PC_FieldSplit  *jac;
23750971522cSBarry Smith 
23760971522cSBarry Smith   PetscFunctionBegin;
2377b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
23782fa5cd67SKarl Rupp 
23790971522cSBarry Smith   jac->bs                 = -1;
23800971522cSBarry Smith   jac->nsplits            = 0;
23813e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2382e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2383c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2384c096484dSStefano Zampini   jac->schurscale         = -1.0;
2385fbe7908bSJed Brown   jac->dm_splits          = PETSC_TRUE;
238651f519a2SBarry Smith 
23870971522cSBarry Smith   pc->data = (void*)jac;
23880971522cSBarry Smith 
23890971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
2390421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
23910971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
2392574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
23930971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
23940971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
23950971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
23960971522cSBarry Smith   pc->ops->applyrichardson = 0;
23970971522cSBarry Smith 
2398bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
2399bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
2400bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
2401bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
2402bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
24036dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);CHKERRQ(ierr);
24040971522cSBarry Smith   PetscFunctionReturn(0);
24050971522cSBarry Smith }
24060971522cSBarry Smith 
24070971522cSBarry Smith 
2408a541d17aSBarry Smith 
2409