xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision de482cd7755a69d9c0e371a0f04bcf31363d03f1)
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 
49a51937d4SCarola Kruse   /* Only used when Golub-Kahan bidiagonalization preconditioning is used */
50a51937d4SCarola Kruse   Mat                       H;                     /* The modified matrix H = A00 + nu*A01*A01'              */
51a51937d4SCarola Kruse   PetscReal                 gkbtol;                /* Stopping tolerance for lower bound estimate            */
52a51937d4SCarola Kruse   PetscInt                  gkbdelay;              /* The delay window for the stopping criterion            */
53a51937d4SCarola Kruse   PetscReal                 gkbnu;                 /* Parameter for augmented Lagrangian H = A + nu*A01*A01' */
54a51937d4SCarola Kruse   PetscInt                  gkbmaxit;              /* Maximum number of iterations for outer loop            */
55*de482cd7SCarola Kruse   PetscBool                 gkbmonitor;            /* Monitor for gkb iterations and the lower bound error   */
56*de482cd7SCarola Kruse   PetscViewer               gkbviewer;             /* Viewer context for gkbmonitor                          */
57a51937d4SCarola Kruse 
5897bbdb24SBarry Smith   PC_FieldSplitLink         head;
596dbb499eSCian Wilson   PetscBool                 isrestrict;             /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
60c1570756SJed Brown   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
614ab8060aSDmitry Karpeev   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
62c84da90fSDmitry Karpeev   PetscBool                 diag_use_amat;          /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
63c84da90fSDmitry Karpeev   PetscBool                 offdiag_use_amat;       /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
647b752e3dSPatrick Sanan   PetscBool                 detect;                 /* Whether to form 2-way split by finding zero diagonal entries */
650971522cSBarry Smith } PC_FieldSplit;
660971522cSBarry Smith 
6716913363SBarry Smith /*
6895452b02SPatrick Sanan     Notes:
6995452b02SPatrick Sanan     there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
7016913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
7116913363SBarry Smith    PC you could change this.
7216913363SBarry Smith */
73084e4875SJed Brown 
74e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
75084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
76084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
77084e4875SJed Brown {
78084e4875SJed Brown   switch (jac->schurpre) {
79084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
80a7476a74SDmitry Karpeev   case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
81e87fab1bSBarry Smith   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
82e74569cdSMatthew G. Knepley   case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
83084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
84084e4875SJed Brown   default:
85084e4875SJed Brown     return jac->schur_user ? jac->schur_user : jac->pmat[1];
86084e4875SJed Brown   }
87084e4875SJed Brown }
88084e4875SJed Brown 
89084e4875SJed Brown 
909804daf3SBarry Smith #include <petscdraw.h>
910971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
920971522cSBarry Smith {
930971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
940971522cSBarry Smith   PetscErrorCode    ierr;
95d9884438SBarry Smith   PetscBool         iascii,isdraw;
960971522cSBarry Smith   PetscInt          i,j;
975a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
980971522cSBarry Smith 
990971522cSBarry Smith   PetscFunctionBegin;
100251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
101d9884438SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1020971522cSBarry Smith   if (iascii) {
1032c9966d7SBarry Smith     if (jac->bs > 0) {
10451f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
1052c9966d7SBarry Smith     } else {
1062c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
1072c9966d7SBarry Smith     }
108f5236f50SJed Brown     if (pc->useAmat) {
109f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
110a3df900dSMatthew G Knepley     }
111c84da90fSDmitry Karpeev     if (jac->diag_use_amat) {
112c84da90fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for diagonal blocks\n");CHKERRQ(ierr);
113c84da90fSDmitry Karpeev     }
114c84da90fSDmitry Karpeev     if (jac->offdiag_use_amat) {
115c84da90fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for off-diagonal blocks\n");CHKERRQ(ierr);
116c84da90fSDmitry Karpeev     }
11769a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
1180971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
1191ab39975SBarry Smith       if (ilink->fields) {
1200971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
12179416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1225a9f2f41SSatish Balay         for (j=0; j<ilink->nfields; j++) {
12379416396SBarry Smith           if (j > 0) {
12479416396SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
12579416396SBarry Smith           }
1265a9f2f41SSatish Balay           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1270971522cSBarry Smith         }
1280971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
12979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1301ab39975SBarry Smith       } else {
1311ab39975SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1321ab39975SBarry Smith       }
1335a9f2f41SSatish Balay       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1345a9f2f41SSatish Balay       ilink = ilink->next;
1350971522cSBarry Smith     }
1362fa5cd67SKarl Rupp   }
1372fa5cd67SKarl Rupp 
1382fa5cd67SKarl Rupp  if (isdraw) {
139d9884438SBarry Smith     PetscDraw draw;
140d9884438SBarry Smith     PetscReal x,y,w,wd;
141d9884438SBarry Smith 
142d9884438SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
143d9884438SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
144d9884438SBarry Smith     w    = 2*PetscMin(1.0 - x,x);
145d9884438SBarry Smith     wd   = w/(jac->nsplits + 1);
146d9884438SBarry Smith     x    = x - wd*(jac->nsplits-1)/2.0;
147d9884438SBarry Smith     for (i=0; i<jac->nsplits; i++) {
148d9884438SBarry Smith       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
149d9884438SBarry Smith       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
150d9884438SBarry Smith       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
151d9884438SBarry Smith       x    += wd;
152d9884438SBarry Smith       ilink = ilink->next;
153d9884438SBarry Smith     }
1540971522cSBarry Smith   }
1550971522cSBarry Smith   PetscFunctionReturn(0);
1560971522cSBarry Smith }
1570971522cSBarry Smith 
1583b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1593b224e63SBarry Smith {
1603b224e63SBarry Smith   PC_FieldSplit              *jac = (PC_FieldSplit*)pc->data;
1613b224e63SBarry Smith   PetscErrorCode             ierr;
1624996c5bdSBarry Smith   PetscBool                  iascii,isdraw;
1633b224e63SBarry Smith   PetscInt                   i,j;
1643b224e63SBarry Smith   PC_FieldSplitLink          ilink = jac->head;
165a9908d51SBarry Smith   MatSchurComplementAinvType atype;
1663b224e63SBarry Smith 
1673b224e63SBarry Smith   PetscFunctionBegin;
168251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1694996c5bdSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1703b224e63SBarry Smith   if (iascii) {
1712c9966d7SBarry Smith     if (jac->bs > 0) {
172c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1732c9966d7SBarry Smith     } else {
1742c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1752c9966d7SBarry Smith     }
176f5236f50SJed Brown     if (pc->useAmat) {
177f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
178a3df900dSMatthew G Knepley     }
1793e8b8b31SMatthew G Knepley     switch (jac->schurpre) {
1803e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
181a9908d51SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);
182a9908d51SBarry Smith       break;
183a7476a74SDmitry Karpeev     case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
184a9908d51SBarry Smith       ierr = MatSchurComplementGetAinvType(jac->schur,&atype);CHKERRQ(ierr);
185d0a9c1a2SBarry 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;
186e87fab1bSBarry Smith     case PC_FIELDSPLIT_SCHUR_PRE_A11:
187a9908d51SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
188a9908d51SBarry Smith       break;
189e74569cdSMatthew G. Knepley     case PC_FIELDSPLIT_SCHUR_PRE_FULL:
190a9908d51SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the exact Schur complement\n");CHKERRQ(ierr);
191a9908d51SBarry Smith       break;
1923e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1933e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1943e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1953e8b8b31SMatthew G Knepley       } else {
196e87fab1bSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
1973e8b8b31SMatthew G Knepley       }
1983e8b8b31SMatthew G Knepley       break;
1993e8b8b31SMatthew G Knepley     default:
20082f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
2013e8b8b31SMatthew G Knepley     }
2023b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
2033b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2043b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
2053b224e63SBarry Smith       if (ilink->fields) {
2063b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
2073b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
2083b224e63SBarry Smith         for (j=0; j<ilink->nfields; j++) {
2093b224e63SBarry Smith           if (j > 0) {
2103b224e63SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
2113b224e63SBarry Smith           }
2123b224e63SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
2133b224e63SBarry Smith         }
2143b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
2153b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
2163b224e63SBarry Smith       } else {
2173b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
2183b224e63SBarry Smith       }
2193b224e63SBarry Smith       ilink = ilink->next;
2203b224e63SBarry Smith     }
221435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");CHKERRQ(ierr);
2223b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
22306de4afeSJed Brown     if (jac->head) {
224443836d0SMatthew G Knepley       ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
22506de4afeSJed Brown     } else  {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
2263b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
22706de4afeSJed Brown     if (jac->head && jac->kspupper != jac->head->ksp) {
228443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr);
229443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
230b2750c55SJed Brown       if (jac->kspupper) {ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);}
231b2750c55SJed Brown       else {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
232443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
233443836d0SMatthew G Knepley     }
234435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
2353b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
23612cae6f2SJed Brown     if (jac->kspschur) {
2373b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
23812cae6f2SJed Brown     } else {
23912cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
24012cae6f2SJed Brown     }
2413b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2423b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
24306de4afeSJed Brown   } else if (isdraw && jac->head) {
2444996c5bdSBarry Smith     PetscDraw draw;
2454996c5bdSBarry Smith     PetscReal x,y,w,wd,h;
2464996c5bdSBarry Smith     PetscInt  cnt = 2;
2474996c5bdSBarry Smith     char      str[32];
2484996c5bdSBarry Smith 
2494996c5bdSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2504996c5bdSBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
251c74581afSBarry Smith     if (jac->kspupper != jac->head->ksp) cnt++;
252c74581afSBarry Smith     w  = 2*PetscMin(1.0 - x,x);
253c74581afSBarry Smith     wd = w/(cnt + 1);
254c74581afSBarry Smith 
2554996c5bdSBarry Smith     ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
25651fa3d41SBarry Smith     ierr = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2574996c5bdSBarry Smith     y   -= h;
2584996c5bdSBarry Smith     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
259e87fab1bSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr);
2603b224e63SBarry Smith     } else {
2614996c5bdSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr);
2624996c5bdSBarry Smith     }
26351fa3d41SBarry Smith     ierr = PetscDrawStringBoxed(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2644996c5bdSBarry Smith     y   -= h;
2654996c5bdSBarry Smith     x    = x - wd*(cnt-1)/2.0;
2664996c5bdSBarry Smith 
2674996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2684996c5bdSBarry Smith     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
2694996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2704996c5bdSBarry Smith     if (jac->kspupper != jac->head->ksp) {
2714996c5bdSBarry Smith       x   += wd;
2724996c5bdSBarry Smith       ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2734996c5bdSBarry Smith       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
2744996c5bdSBarry Smith       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2754996c5bdSBarry Smith     }
2764996c5bdSBarry Smith     x   += wd;
2774996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2784996c5bdSBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
2794996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2803b224e63SBarry Smith   }
2813b224e63SBarry Smith   PetscFunctionReturn(0);
2823b224e63SBarry Smith }
2833b224e63SBarry Smith 
284a51937d4SCarola Kruse static PetscErrorCode PCView_FieldSplit_GKB(PC pc,PetscViewer viewer)
285a51937d4SCarola Kruse {
286a51937d4SCarola Kruse   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
287a51937d4SCarola Kruse   PetscErrorCode    ierr;
288a51937d4SCarola Kruse   PetscBool         iascii,isdraw;
289a51937d4SCarola Kruse   PetscInt          i,j;
290a51937d4SCarola Kruse   PC_FieldSplitLink ilink = jac->head;
291a51937d4SCarola Kruse 
292a51937d4SCarola Kruse   PetscFunctionBegin;
293a51937d4SCarola Kruse   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
294a51937d4SCarola Kruse   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
295a51937d4SCarola Kruse   if (iascii) {
296a51937d4SCarola Kruse     if (jac->bs > 0) {
297a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
298a51937d4SCarola Kruse     } else {
299a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
300a51937d4SCarola Kruse     }
301a51937d4SCarola Kruse     if (pc->useAmat) {
302a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
303a51937d4SCarola Kruse     }
304a51937d4SCarola Kruse     if (jac->diag_use_amat) {
305a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for diagonal blocks\n");CHKERRQ(ierr);
306a51937d4SCarola Kruse     }
307a51937d4SCarola Kruse     if (jac->offdiag_use_amat) {
308a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for off-diagonal blocks\n");CHKERRQ(ierr);
309a51937d4SCarola Kruse     }
310a51937d4SCarola Kruse 
311a51937d4SCarola Kruse     ierr = PetscViewerASCIIPrintf(viewer,"  Stopping tolerance=%.1e, delay in error estimate=%D, maximum iterations=%D\n",jac->gkbtol,jac->gkbdelay,jac->gkbmaxit);CHKERRQ(ierr);
312a51937d4SCarola Kruse     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for H = A00 + nu*A01*A01' matrix:\n");CHKERRQ(ierr);
313a51937d4SCarola Kruse     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
314a51937d4SCarola Kruse 
315a51937d4SCarola Kruse     if (ilink->fields) {
316a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",0);CHKERRQ(ierr);
317a51937d4SCarola Kruse       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
318a51937d4SCarola Kruse       for (j=0; j<ilink->nfields; j++) {
319a51937d4SCarola Kruse         if (j > 0) {
320a51937d4SCarola Kruse           ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
321a51937d4SCarola Kruse         }
322a51937d4SCarola Kruse         ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
323a51937d4SCarola Kruse       }
324a51937d4SCarola Kruse       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
325a51937d4SCarola Kruse       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
326a51937d4SCarola Kruse     } else {
327a51937d4SCarola Kruse         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",0);CHKERRQ(ierr);
328a51937d4SCarola Kruse     }
329a51937d4SCarola Kruse     ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
330a51937d4SCarola Kruse 
331a51937d4SCarola Kruse     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
332a51937d4SCarola Kruse   }
333a51937d4SCarola Kruse 
334a51937d4SCarola Kruse  if (isdraw) {
335a51937d4SCarola Kruse     PetscDraw draw;
336a51937d4SCarola Kruse     PetscReal x,y,w,wd;
337a51937d4SCarola Kruse 
338a51937d4SCarola Kruse     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
339a51937d4SCarola Kruse     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
340a51937d4SCarola Kruse     w    = 2*PetscMin(1.0 - x,x);
341a51937d4SCarola Kruse     wd   = w/(jac->nsplits + 1);
342a51937d4SCarola Kruse     x    = x - wd*(jac->nsplits-1)/2.0;
343a51937d4SCarola Kruse     for (i=0; i<jac->nsplits; i++) {
344a51937d4SCarola Kruse       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
345a51937d4SCarola Kruse       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
346a51937d4SCarola Kruse       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
347a51937d4SCarola Kruse       x    += wd;
348a51937d4SCarola Kruse       ilink = ilink->next;
349a51937d4SCarola Kruse     }
350a51937d4SCarola Kruse   }
351a51937d4SCarola Kruse   PetscFunctionReturn(0);
352a51937d4SCarola Kruse }
353a51937d4SCarola Kruse 
354a51937d4SCarola Kruse 
3556c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
3566c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
3576c924f48SJed Brown {
3586c924f48SJed Brown   PetscErrorCode ierr;
3596c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
3605d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
3615d4c12cdSJungho Lee   PetscBool      flg,flg_col;
3625d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
3636c924f48SJed Brown 
3646c924f48SJed Brown   PetscFunctionBegin;
365785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields);CHKERRQ(ierr);
366785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields_col);CHKERRQ(ierr);
3676c924f48SJed Brown   for (i=0,flg=PETSC_TRUE;; i++) {
3688caf3d72SBarry Smith     ierr        = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
3698caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
3708caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
3716c924f48SJed Brown     nfields     = jac->bs;
37229499fbbSJungho Lee     nfields_col = jac->bs;
373c5929fdfSBarry Smith     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
374c5929fdfSBarry Smith     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
3756c924f48SJed Brown     if (!flg) break;
3765d4c12cdSJungho Lee     else if (flg && !flg_col) {
3775d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
3785d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
3792fa5cd67SKarl Rupp     } else {
3805d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
3815d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
3825d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
3835d4c12cdSJungho Lee     }
3846c924f48SJed Brown   }
3856c924f48SJed Brown   if (i > 0) {
3866c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
3876c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
3886c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
3896c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
3906c924f48SJed Brown   }
3916c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
3925d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
3936c924f48SJed Brown   PetscFunctionReturn(0);
3946c924f48SJed Brown }
3956c924f48SJed Brown 
39669a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
3970971522cSBarry Smith {
3980971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3990971522cSBarry Smith   PetscErrorCode    ierr;
4005a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
4017b752e3dSPatrick Sanan   PetscBool         fieldsplit_default = PETSC_FALSE,coupling = PETSC_FALSE;
4026c924f48SJed Brown   PetscInt          i;
4030971522cSBarry Smith 
4040971522cSBarry Smith   PetscFunctionBegin;
4057287d2a3SDmitry Karpeev   /*
406f5f0d762SBarry Smith    Kinda messy, but at least this now uses DMCreateFieldDecomposition().
4077287d2a3SDmitry Karpeev    Should probably be rewritten.
4087287d2a3SDmitry Karpeev    */
409f5f0d762SBarry Smith   if (!ilink) {
410c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_coupling",&coupling,NULL);CHKERRQ(ierr);
4117b752e3dSPatrick Sanan     if (pc->dm && jac->dm_splits && !jac->detect && !coupling) {
412bafc1b83SMatthew G Knepley       PetscInt  numFields, f, i, j;
4130784a22cSJed Brown       char      **fieldNames;
4147b62db95SJungho Lee       IS        *fields;
415e7c4fc90SDmitry Karpeev       DM        *dms;
416bafc1b83SMatthew G Knepley       DM        subdm[128];
417bafc1b83SMatthew G Knepley       PetscBool flg;
418bafc1b83SMatthew G Knepley 
41916621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
420bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
421bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE;; i++) {
422bafc1b83SMatthew G Knepley         PetscInt ifields[128];
423bafc1b83SMatthew G Knepley         IS       compField;
424bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
425bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
426bafc1b83SMatthew G Knepley 
4278caf3d72SBarry Smith         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
428c5929fdfSBarry Smith         ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
429bafc1b83SMatthew G Knepley         if (!flg) break;
43082f516ccSBarry Smith         if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
431bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
432bafc1b83SMatthew G Knepley         if (nfields == 1) {
433bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
434bafc1b83SMatthew G Knepley         } else {
4358caf3d72SBarry Smith           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
436bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
4377287d2a3SDmitry Karpeev         }
438bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
439bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
440bafc1b83SMatthew G Knepley           f    = ifields[j];
4417b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
4427b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
4437b62db95SJungho Lee         }
444bafc1b83SMatthew G Knepley       }
445bafc1b83SMatthew G Knepley       if (i == 0) {
446bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
447bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
448bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
449bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
450bafc1b83SMatthew G Knepley         }
451bafc1b83SMatthew G Knepley       } else {
452d724dfffSBarry Smith         for (j=0; j<numFields; j++) {
453d724dfffSBarry Smith           ierr = DMDestroy(dms+j);CHKERRQ(ierr);
454d724dfffSBarry Smith         }
455d724dfffSBarry Smith         ierr = PetscFree(dms);CHKERRQ(ierr);
456785e854fSJed Brown         ierr = PetscMalloc1(i, &dms);CHKERRQ(ierr);
4572fa5cd67SKarl Rupp         for (j = 0; j < i; ++j) dms[j] = subdm[j];
458bafc1b83SMatthew G Knepley       }
4597b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
4607b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
461e7c4fc90SDmitry Karpeev       if (dms) {
4628b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
463bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
4647287d2a3SDmitry Karpeev           const char *prefix;
4657287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr);
4667287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr);
4677b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
4687b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
4697287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr);
470e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
4712fa5ba8aSJed Brown         }
4727b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
4738b8307b2SJed Brown       }
47466ffff09SJed Brown     } else {
475521d7252SBarry Smith       if (jac->bs <= 0) {
476704ba839SBarry Smith         if (pc->pmat) {
477521d7252SBarry Smith           ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
4782fa5cd67SKarl Rupp         } else jac->bs = 1;
479521d7252SBarry Smith       }
480d32f9abdSBarry Smith 
4817b752e3dSPatrick Sanan       if (jac->detect) {
4826ce1633cSBarry Smith         IS       zerodiags,rest;
4836ce1633cSBarry Smith         PetscInt nmin,nmax;
4846ce1633cSBarry Smith 
4856ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4866ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
4876ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
4886ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
4896ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
490fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
491fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4923a062f41SBarry Smith       } else if (coupling) {
4933a062f41SBarry Smith         IS       coupling,rest;
4943a062f41SBarry Smith         PetscInt nmin,nmax;
4953a062f41SBarry Smith 
4963a062f41SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4973a062f41SBarry Smith         ierr = MatFindOffBlockDiagonalEntries(pc->mat,&coupling);CHKERRQ(ierr);
4986bea0878SBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);CHKERRQ(ierr);
499e52d2c62SBarry Smith         ierr = ISSetIdentity(rest);CHKERRQ(ierr);
500d020c841SBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
501d020c841SBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",coupling);CHKERRQ(ierr);
5023a062f41SBarry Smith         ierr = ISDestroy(&coupling);CHKERRQ(ierr);
5033a062f41SBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
5046ce1633cSBarry Smith       } else {
505c5929fdfSBarry Smith         ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
5067287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
507d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
508d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
5096c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
5106c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
511d32f9abdSBarry Smith         }
5126dbb499eSCian Wilson         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
513d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
514db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
5156c924f48SJed Brown             char splitname[8];
5168caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
5175d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
51879416396SBarry Smith           }
5195d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
520521d7252SBarry Smith         }
52166ffff09SJed Brown       }
5226ce1633cSBarry Smith     }
523edf189efSBarry Smith   } else if (jac->nsplits == 1) {
524edf189efSBarry Smith     if (ilink->is) {
525edf189efSBarry Smith       IS       is2;
526edf189efSBarry Smith       PetscInt nmin,nmax;
527edf189efSBarry Smith 
528edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
529edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
530db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
531fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
53282f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
533edf189efSBarry Smith   }
534d0af7cd3SBarry Smith 
53582f516ccSBarry Smith   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
53669a612a9SBarry Smith   PetscFunctionReturn(0);
53769a612a9SBarry Smith }
53869a612a9SBarry Smith 
539a51937d4SCarola Kruse static PetscErrorCode MatGolubKahanComputeExplicitOperator(Mat A,Mat B,Mat C,Mat *H,PetscReal gkbnu)
540a51937d4SCarola Kruse {
541a51937d4SCarola Kruse   PetscErrorCode    ierr;
542a51937d4SCarola Kruse   Mat               BT,T;
543*de482cd7SCarola Kruse   PetscReal         nrmT,nrmB;
544a51937d4SCarola Kruse 
545a51937d4SCarola Kruse   PetscFunctionBegin;
546a51937d4SCarola Kruse   ierr = MatHermitianTranspose(C,MAT_INITIAL_MATRIX,&T);CHKERRQ(ierr);            /* Test if augmented matrix is symmetric */
547a51937d4SCarola Kruse   ierr = MatAXPY(T,-1.0,B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
548*de482cd7SCarola Kruse   ierr = MatNorm(T,NORM_1,&nrmT);CHKERRQ(ierr);
549*de482cd7SCarola Kruse   ierr = MatNorm(B,NORM_1,&nrmB);CHKERRQ(ierr);
550*de482cd7SCarola Kruse   if (nrmB > 0) {
551*de482cd7SCarola Kruse     if (nrmT/nrmB >= PETSC_SMALL) {
552*de482cd7SCarola Kruse       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix is not symmetric/hermitian, GKB is not applicable.");
553*de482cd7SCarola Kruse     }
554a51937d4SCarola Kruse   }
555a51937d4SCarola Kruse   /* Compute augmented Lagrangian matrix H = A00 + nu*A01*A01'. This corresponds to */
556a51937d4SCarola Kruse   /* setting N := 1/nu*I in [Ar13].                                                 */
557a51937d4SCarola Kruse   ierr = MatHermitianTranspose(B,MAT_INITIAL_MATRIX,&BT);CHKERRQ(ierr);
558a51937d4SCarola Kruse   ierr = MatMatMult(B,BT,MAT_INITIAL_MATRIX,PETSC_DEFAULT,H);CHKERRQ(ierr);       /* H = A01*A01'          */
559a51937d4SCarola Kruse   ierr = MatAYPX(*H,gkbnu,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);             /* H = A00 + nu*A01*A01' */
560a51937d4SCarola Kruse 
561a51937d4SCarola Kruse   ierr = MatDestroy(&BT);CHKERRQ(ierr);
562a51937d4SCarola Kruse   ierr = MatDestroy(&T);CHKERRQ(ierr);
563a51937d4SCarola Kruse   PetscFunctionReturn(0);
564a51937d4SCarola Kruse }
565a51937d4SCarola Kruse 
5662d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions,const char pre[], const char name[],const char *value[],PetscBool *flg);
567514bf10dSMatthew G Knepley 
56869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
56969a612a9SBarry Smith {
57069a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
57169a612a9SBarry Smith   PetscErrorCode    ierr;
5725a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
5732c9966d7SBarry Smith   PetscInt          i,nsplit;
5745f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
57569a612a9SBarry Smith 
57669a612a9SBarry Smith   PetscFunctionBegin;
57769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
57897bbdb24SBarry Smith   nsplit = jac->nsplits;
5795a9f2f41SSatish Balay   ilink  = jac->head;
58097bbdb24SBarry Smith 
58197bbdb24SBarry Smith   /* get the matrices for each split */
582704ba839SBarry Smith   if (!jac->issetup) {
5831b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
58497bbdb24SBarry Smith 
585704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
586704ba839SBarry Smith 
5875d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
5882c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
5892c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
5902c9966d7SBarry Smith     }
59151f519a2SBarry Smith     bs     = jac->bs;
59297bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
5931b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
5941b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
5951b9fc7fcSBarry Smith       if (jac->defaultsplit) {
596ce94432eSBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
5975f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
598704ba839SBarry Smith       } else if (!ilink->is) {
599ccb205f8SBarry Smith         if (ilink->nfields > 1) {
6005f4ab4e1SJungho Lee           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
601785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&ii);CHKERRQ(ierr);
602785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&jj);CHKERRQ(ierr);
6031b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
6041b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
6051b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
6065f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
60797bbdb24SBarry Smith             }
60897bbdb24SBarry Smith           }
609ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
610ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
61190e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is, nfields);CHKERRQ(ierr);
61290e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is_col, nfields);CHKERRQ(ierr);
613ccb205f8SBarry Smith         } else {
614ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
615ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
616ccb205f8SBarry Smith         }
6173e197d65SBarry Smith       }
618704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
6195f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
6205f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
621704ba839SBarry Smith       ilink = ilink->next;
6221b9fc7fcSBarry Smith     }
6231b9fc7fcSBarry Smith   }
6241b9fc7fcSBarry Smith 
625704ba839SBarry Smith   ilink = jac->head;
62697bbdb24SBarry Smith   if (!jac->pmat) {
627bdddcaaaSMatthew G Knepley     Vec xtmp;
628bdddcaaaSMatthew G Knepley 
6292a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
630785e854fSJed Brown     ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr);
631dcca6d9dSJed Brown     ierr = PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);CHKERRQ(ierr);
632cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
633bdddcaaaSMatthew G Knepley       MatNullSpace sp;
634bdddcaaaSMatthew G Knepley 
635a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
636a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
637a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
638a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
639a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
640a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
6412fa5cd67SKarl Rupp 
642a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
643a3df900dSMatthew G Knepley         }
644a3df900dSMatthew G Knepley       } else {
6453a062f41SBarry Smith         const char *prefix;
6467dae84e0SHong Zhang         ierr = MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
6473a062f41SBarry Smith         ierr = KSPGetOptionsPrefix(ilink->ksp,&prefix);CHKERRQ(ierr);
6483a062f41SBarry Smith         ierr = MatSetOptionsPrefix(jac->pmat[i],prefix);CHKERRQ(ierr);
6493a062f41SBarry Smith         ierr = MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");CHKERRQ(ierr);
650a3df900dSMatthew G Knepley       }
651bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
6522a7a6963SBarry Smith       ierr     = MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
6530298fd71SBarry Smith       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
654bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
65535928de7SBarry Smith       ierr = VecScatterCreateWithData(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
656ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
657ed1f0337SMatthew G Knepley       if (sp) {
658ed1f0337SMatthew G Knepley         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
659ed1f0337SMatthew G Knepley       }
660704ba839SBarry Smith       ilink = ilink->next;
661cf502942SBarry Smith     }
662bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
66397bbdb24SBarry Smith   } else {
664ef7efd37SHong Zhang     MatReuse scall;
665ef7efd37SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
666ef7efd37SHong Zhang       for (i=0; i<nsplit; i++) {
667ef7efd37SHong Zhang         ierr = MatDestroy(&jac->pmat[i]);CHKERRQ(ierr);
668ef7efd37SHong Zhang       }
669ef7efd37SHong Zhang       scall = MAT_INITIAL_MATRIX;
670ef7efd37SHong Zhang     } else scall = MAT_REUSE_MATRIX;
671ef7efd37SHong Zhang 
672cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
673a3df900dSMatthew G Knepley       Mat pmat;
674a3df900dSMatthew G Knepley 
675a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
676a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
677a3df900dSMatthew G Knepley       if (!pmat) {
678ef7efd37SHong Zhang         ierr = MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,scall,&jac->pmat[i]);CHKERRQ(ierr);
679a3df900dSMatthew G Knepley       }
680704ba839SBarry Smith       ilink = ilink->next;
681cf502942SBarry Smith     }
68297bbdb24SBarry Smith   }
6834e39094bSDmitry Karpeev   if (jac->diag_use_amat) {
684519d70e2SJed Brown     ilink = jac->head;
685519d70e2SJed Brown     if (!jac->mat) {
686785e854fSJed Brown       ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr);
687519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
6887dae84e0SHong Zhang         ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
689519d70e2SJed Brown         ilink = ilink->next;
690519d70e2SJed Brown       }
691519d70e2SJed Brown     } else {
692ef7efd37SHong Zhang       MatReuse scall;
693ef7efd37SHong Zhang       if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
694519d70e2SJed Brown         for (i=0; i<nsplit; i++) {
695ef7efd37SHong Zhang           ierr = MatDestroy(&jac->mat[i]);CHKERRQ(ierr);
696ef7efd37SHong Zhang         }
697ef7efd37SHong Zhang         scall = MAT_INITIAL_MATRIX;
698ef7efd37SHong Zhang       } else scall = MAT_REUSE_MATRIX;
699ef7efd37SHong Zhang 
700ef7efd37SHong Zhang       for (i=0; i<nsplit; i++) {
701ef7efd37SHong Zhang         if (jac->mat[i]) {ierr = MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,scall,&jac->mat[i]);CHKERRQ(ierr);}
702519d70e2SJed Brown         ilink = ilink->next;
703519d70e2SJed Brown       }
704519d70e2SJed Brown     }
705519d70e2SJed Brown   } else {
706519d70e2SJed Brown     jac->mat = jac->pmat;
707519d70e2SJed Brown   }
70897bbdb24SBarry Smith 
70953935eafSBarry Smith   /* Check for null space attached to IS */
71053935eafSBarry Smith   ilink = jac->head;
71153935eafSBarry Smith   for (i=0; i<nsplit; i++) {
71253935eafSBarry Smith     MatNullSpace sp;
71353935eafSBarry Smith 
71453935eafSBarry Smith     ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
71553935eafSBarry Smith     if (sp) {
71653935eafSBarry Smith       ierr = MatSetNullSpace(jac->mat[i], sp);CHKERRQ(ierr);
71753935eafSBarry Smith     }
71853935eafSBarry Smith     ilink = ilink->next;
71953935eafSBarry Smith   }
72053935eafSBarry Smith 
721a51937d4SCarola Kruse   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR && jac->type != PC_COMPOSITE_GKB) {
72268dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
7234e39094bSDmitry Karpeev     /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
72468dd23aaSBarry Smith     ilink = jac->head;
725e52d2c62SBarry Smith     if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
726e52d2c62SBarry 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 */
727e52d2c62SBarry Smith       if (!jac->Afield) {
728e52d2c62SBarry Smith         ierr = PetscCalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
72980c96bb1SFande Kong         if (jac->offdiag_use_amat) {
7307dae84e0SHong Zhang           ierr  = MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
731e52d2c62SBarry Smith         } else {
7327dae84e0SHong Zhang           ierr  = MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
73380c96bb1SFande Kong         }
73480c96bb1SFande Kong       } else {
735ef7efd37SHong Zhang         MatReuse scall;
736ef7efd37SHong Zhang         if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
737ef7efd37SHong Zhang           for (i=0; i<nsplit; i++) {
738ef7efd37SHong Zhang             ierr = MatDestroy(&jac->Afield[1]);CHKERRQ(ierr);
739ef7efd37SHong Zhang           }
740ef7efd37SHong Zhang           scall = MAT_INITIAL_MATRIX;
741ef7efd37SHong Zhang         } else scall = MAT_REUSE_MATRIX;
742ef7efd37SHong Zhang 
74380c96bb1SFande Kong         if (jac->offdiag_use_amat) {
744ef7efd37SHong Zhang           ierr  = MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,scall,&jac->Afield[1]);CHKERRQ(ierr);
74580c96bb1SFande Kong         } else {
746ef7efd37SHong Zhang           ierr  = MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,scall,&jac->Afield[1]);CHKERRQ(ierr);
74780c96bb1SFande Kong         }
748e52d2c62SBarry Smith       }
749e52d2c62SBarry Smith     } else {
75068dd23aaSBarry Smith       if (!jac->Afield) {
751785e854fSJed Brown         ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
75268dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
75380c96bb1SFande Kong           if (jac->offdiag_use_amat) {
7547dae84e0SHong Zhang             ierr  = MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
75580c96bb1SFande Kong           } else {
7567dae84e0SHong Zhang             ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
75780c96bb1SFande Kong           }
75868dd23aaSBarry Smith           ilink = ilink->next;
75968dd23aaSBarry Smith         }
76068dd23aaSBarry Smith       } else {
761ef7efd37SHong Zhang         MatReuse scall;
762ef7efd37SHong Zhang         if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
763ef7efd37SHong Zhang           for (i=0; i<nsplit; i++) {
764ef7efd37SHong Zhang             ierr = MatDestroy(&jac->Afield[i]);CHKERRQ(ierr);
765ef7efd37SHong Zhang           }
766ef7efd37SHong Zhang           scall = MAT_INITIAL_MATRIX;
767ef7efd37SHong Zhang         } else scall = MAT_REUSE_MATRIX;
768ef7efd37SHong Zhang 
76968dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
77080c96bb1SFande Kong           if (jac->offdiag_use_amat) {
771ef7efd37SHong Zhang             ierr  = MatCreateSubMatrix(pc->mat,ilink->is,NULL,scall,&jac->Afield[i]);CHKERRQ(ierr);
77280c96bb1SFande Kong           } else {
773ef7efd37SHong Zhang             ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,NULL,scall,&jac->Afield[i]);CHKERRQ(ierr);
77480c96bb1SFande Kong           }
77568dd23aaSBarry Smith           ilink = ilink->next;
77668dd23aaSBarry Smith         }
77768dd23aaSBarry Smith       }
77868dd23aaSBarry Smith     }
779e52d2c62SBarry Smith   }
78068dd23aaSBarry Smith 
7813b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
7823b224e63SBarry Smith     IS          ccis;
783c096484dSStefano Zampini     PetscBool   isspd;
7844aa3045dSJed Brown     PetscInt    rstart,rend;
785093c86ffSJed Brown     char        lscname[256];
786093c86ffSJed Brown     PetscObject LSC_L;
787ce94432eSBarry Smith 
788ce94432eSBarry Smith     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
78968dd23aaSBarry Smith 
790c096484dSStefano Zampini     /* If pc->mat is SPD, don't scale by -1 the Schur complement */
791c096484dSStefano Zampini     if (jac->schurscale == (PetscScalar)-1.0) {
792c096484dSStefano Zampini       ierr = MatGetOption(pc->pmat,MAT_SPD,&isspd);CHKERRQ(ierr);
793c096484dSStefano Zampini       jac->schurscale = (isspd == PETSC_TRUE) ? 1.0 : -1.0;
794c096484dSStefano Zampini     }
795c096484dSStefano Zampini 
796e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
797e6cab6aaSJed Brown     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
798e6cab6aaSJed Brown 
7993b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
8003b224e63SBarry Smith     if (jac->schur) {
8010298fd71SBarry Smith       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
802443836d0SMatthew G Knepley 
803fb3147dbSMatthew G Knepley       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
8043b224e63SBarry Smith       ilink = jac->head;
80549bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
8064e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8077dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
808c84da90fSDmitry Karpeev       } else {
8097dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
810c84da90fSDmitry Karpeev       }
811fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
8123b224e63SBarry Smith       ilink = ilink->next;
81349bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
8144e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8157dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
816c84da90fSDmitry Karpeev       } else {
8177dae84e0SHong Zhang 	ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
818c84da90fSDmitry Karpeev       }
819fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
820aa6c7ce3SBarry Smith       ierr  = MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
821a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
822a7476a74SDmitry Karpeev 	ierr = MatDestroy(&jac->schurp);CHKERRQ(ierr);
823a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
824a7476a74SDmitry Karpeev       }
825443836d0SMatthew G Knepley       if (kspA != kspInner) {
82623ee1639SBarry Smith         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
827443836d0SMatthew G Knepley       }
828443836d0SMatthew G Knepley       if (kspUpper != kspA) {
82923ee1639SBarry Smith         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
830443836d0SMatthew G Knepley       }
83123ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
8323b224e63SBarry Smith     } else {
833bafc1b83SMatthew G Knepley       const char   *Dprefix;
834470b340bSDmitry Karpeev       char         schurprefix[256], schurmatprefix[256];
835514bf10dSMatthew G Knepley       char         schurtestoption[256];
836bdddcaaaSMatthew G Knepley       MatNullSpace sp;
837514bf10dSMatthew G Knepley       PetscBool    flg;
838686bed4dSStefano Zampini       KSP          kspt;
8393b224e63SBarry Smith 
840a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
8413b224e63SBarry Smith       ilink = jac->head;
84249bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
8434e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8447dae84e0SHong Zhang 	ierr = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
845c84da90fSDmitry Karpeev       } else {
8467dae84e0SHong Zhang 	ierr = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
847c84da90fSDmitry Karpeev       }
848fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
8493b224e63SBarry Smith       ilink = ilink->next;
85049bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
8514e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8527dae84e0SHong Zhang 	ierr = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
853c84da90fSDmitry Karpeev       } else {
8547dae84e0SHong Zhang 	ierr = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
855c84da90fSDmitry Karpeev       }
856fcfd50ebSBarry Smith       ierr = ISDestroy(&ccis);CHKERRQ(ierr);
85720252d06SBarry Smith 
858f5236f50SJed Brown       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
85920252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
86020252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
861bee83525SDmitry Karpeev       ierr = MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
862470b340bSDmitry Karpeev       ierr = PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
863470b340bSDmitry Karpeev       ierr = MatSetOptionsPrefix(jac->schur,schurmatprefix);CHKERRQ(ierr);
864686bed4dSStefano Zampini       ierr = MatSchurComplementGetKSP(jac->schur,&kspt);CHKERRQ(ierr);
865686bed4dSStefano Zampini       ierr = KSPSetOptionsPrefix(kspt,schurmatprefix);CHKERRQ(ierr);
866686bed4dSStefano Zampini 
867686bed4dSStefano Zampini       /* Note: this is not true in general */
868895c21f2SBarry Smith       ierr = MatGetNullSpace(jac->mat[1], &sp);CHKERRQ(ierr);
86920252d06SBarry Smith       if (sp) {
87020252d06SBarry Smith         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
87120252d06SBarry Smith       }
87220252d06SBarry Smith 
87320252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
874c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
875514bf10dSMatthew G Knepley       if (flg) {
876514bf10dSMatthew G Knepley         DM  dmInner;
87721635b76SJed Brown         KSP kspInner;
878686bed4dSStefano Zampini         PC  pcInner;
87921635b76SJed Brown 
88021635b76SJed Brown         ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
8812cc093acSMatthew G. Knepley         ierr = KSPReset(kspInner);CHKERRQ(ierr);
8822cc093acSMatthew G. Knepley         ierr = KSPSetOperators(kspInner,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
88321635b76SJed Brown         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
88421635b76SJed Brown         /* Indent this deeper to emphasize the "inner" nature of this solver. */
88521635b76SJed Brown         ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr);
8862cc093acSMatthew G. Knepley         ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner->pc, (PetscObject) pc, 2);CHKERRQ(ierr);
88721635b76SJed Brown         ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr);
888514bf10dSMatthew G Knepley 
889514bf10dSMatthew G Knepley         /* Set DM for new solver */
890bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
89121635b76SJed Brown         ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr);
89221635b76SJed Brown         ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr);
893686bed4dSStefano Zampini 
894686bed4dSStefano Zampini         /* Defaults to PCKSP as preconditioner */
895686bed4dSStefano Zampini         ierr = KSPGetPC(kspInner, &pcInner);CHKERRQ(ierr);
896686bed4dSStefano Zampini         ierr = PCSetType(pcInner, PCKSP);CHKERRQ(ierr);
897686bed4dSStefano Zampini         ierr = PCKSPSetKSP(pcInner, jac->head->ksp);CHKERRQ(ierr);
898514bf10dSMatthew G Knepley       } else {
89921635b76SJed Brown          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
90021635b76SJed Brown           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
90121635b76SJed Brown           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
90221635b76SJed 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
90321635b76SJed Brown           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
90421635b76SJed Brown           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
90521635b76SJed Brown         ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr);
906514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr);
907bafc1b83SMatthew G Knepley       }
90823ee1639SBarry Smith       ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
9095a9f2f41SSatish Balay       ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
91097bbdb24SBarry Smith       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
91197bbdb24SBarry Smith 
912686bed4dSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)jac->schur, MATSCHURCOMPLEMENT, &flg);CHKERRQ(ierr);
913686bed4dSStefano Zampini       if (flg) { /* Need to do this otherwise PCSetUp_KSP will overwrite the amat of jac->head->ksp */
914686bed4dSStefano Zampini         KSP kspInner;
915686bed4dSStefano Zampini         PC  pcInner;
916686bed4dSStefano Zampini 
917686bed4dSStefano Zampini         ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
918686bed4dSStefano Zampini         ierr = KSPGetPC(kspInner, &pcInner);CHKERRQ(ierr);
919686bed4dSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)pcInner, PCKSP, &flg);CHKERRQ(ierr);
920686bed4dSStefano Zampini         if (flg) {
921686bed4dSStefano Zampini           KSP ksp;
922686bed4dSStefano Zampini 
923686bed4dSStefano Zampini           ierr = PCKSPGetKSP(pcInner, &ksp);CHKERRQ(ierr);
924686bed4dSStefano Zampini           if (ksp == jac->head->ksp) {
925686bed4dSStefano Zampini             ierr = PCSetUseAmat(pcInner, PETSC_TRUE);CHKERRQ(ierr);
926686bed4dSStefano Zampini           }
927686bed4dSStefano Zampini         }
928686bed4dSStefano Zampini       }
929443836d0SMatthew G Knepley       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
930c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
931443836d0SMatthew G Knepley       if (flg) {
932443836d0SMatthew G Knepley         DM dmInner;
933443836d0SMatthew G Knepley 
934443836d0SMatthew G Knepley         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
93582f516ccSBarry Smith         ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr);
936422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);CHKERRQ(ierr);
937443836d0SMatthew G Knepley         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
9382cc093acSMatthew G. Knepley         ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspupper, (PetscObject) pc, 1);CHKERRQ(ierr);
9392cc093acSMatthew G. Knepley         ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspupper->pc, (PetscObject) pc, 1);CHKERRQ(ierr);
940443836d0SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
941443836d0SMatthew G Knepley         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
942443836d0SMatthew G Knepley         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
943443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
94423ee1639SBarry Smith         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
945443836d0SMatthew G Knepley         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
946443836d0SMatthew G Knepley       } else {
947443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
948443836d0SMatthew G Knepley         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
949443836d0SMatthew G Knepley       }
950443836d0SMatthew G Knepley 
951a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
952a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
953a7476a74SDmitry Karpeev       }
954ce94432eSBarry Smith       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
955422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);CHKERRQ(ierr);
95697bbdb24SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
95720252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
95897bbdb24SBarry Smith       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
9597233a360SDmitry Karpeev         PC pcschur;
9607233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
9617233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
96297bbdb24SBarry Smith         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
963e74569cdSMatthew G. Knepley       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
964e74569cdSMatthew G. Knepley         ierr = MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);CHKERRQ(ierr);
96597bbdb24SBarry Smith       }
96623ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
96797bbdb24SBarry Smith       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
96897bbdb24SBarry Smith       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
969c096484dSStefano Zampini       /* propagate DM */
970b20b4189SMatthew G. Knepley       {
971b20b4189SMatthew G. Knepley         DM sdm;
972b20b4189SMatthew G. Knepley         ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr);
973b20b4189SMatthew G. Knepley         if (sdm) {
974b20b4189SMatthew G. Knepley           ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr);
975b20b4189SMatthew G. Knepley           ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr);
976b20b4189SMatthew G. Knepley         }
977b20b4189SMatthew G. Knepley       }
97897bbdb24SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
97997bbdb24SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
98097bbdb24SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
98197bbdb24SBarry Smith     }
982686bed4dSStefano Zampini     ierr = MatAssemblyBegin(jac->schur,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
983686bed4dSStefano Zampini     ierr = MatAssemblyEnd(jac->schur,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98497bbdb24SBarry Smith 
9855a9f2f41SSatish Balay     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
9868caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
987519d70e2SJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
9883b224e63SBarry Smith     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
989c1570756SJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
9908caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
99197bbdb24SBarry Smith     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
9925a9f2f41SSatish Balay     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
9930971522cSBarry Smith     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
994a51937d4SCarola Kruse   } else if (jac->type == PC_COMPOSITE_GKB) {
995a51937d4SCarola Kruse     IS          ccis;
996a51937d4SCarola Kruse     PetscInt    rstart,rend;
997a51937d4SCarola Kruse 
998a51937d4SCarola Kruse     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use GKB preconditioner you must have exactly 2 fields");
999a51937d4SCarola Kruse 
1000a51937d4SCarola Kruse     ilink = jac->head;
1001a51937d4SCarola Kruse 
1002a51937d4SCarola Kruse     /* When extracting off-diagonal submatrices, we take complements from this range */
1003a51937d4SCarola Kruse     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
1004a51937d4SCarola Kruse 
1005a51937d4SCarola Kruse     ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
1006a51937d4SCarola Kruse     if (jac->offdiag_use_amat) {
1007a51937d4SCarola Kruse      ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
1008a51937d4SCarola Kruse     } else {
1009a51937d4SCarola Kruse       ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
1010a51937d4SCarola Kruse     }
1011a51937d4SCarola Kruse     ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
1012a51937d4SCarola Kruse     ilink = ilink->next;
1013a51937d4SCarola Kruse     ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
1014a51937d4SCarola Kruse     if (jac->offdiag_use_amat) {
1015a51937d4SCarola Kruse       ierr  = MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
1016a51937d4SCarola Kruse     } else {
1017a51937d4SCarola Kruse 	    ierr  = MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
1018a51937d4SCarola Kruse     }
1019a51937d4SCarola Kruse     ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
1020a51937d4SCarola Kruse     ierr  = MatGolubKahanComputeExplicitOperator(jac->mat[0],jac->B,jac->C,&jac->H,jac->gkbnu);CHKERRQ(ierr);
1021a51937d4SCarola Kruse     ilink = jac->head;
1022a51937d4SCarola Kruse     ierr  = KSPSetOperators(ilink->ksp,jac->H,jac->H);CHKERRQ(ierr);
1023a51937d4SCarola Kruse     if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
1024*de482cd7SCarola Kruse     if (jac->gkbmonitor) {
1025*de482cd7SCarola Kruse       PetscInt  tablevel;
1026*de482cd7SCarola Kruse       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&jac->gkbviewer);CHKERRQ(ierr);
1027*de482cd7SCarola Kruse       ierr = PetscViewerSetType(jac->gkbviewer,PETSCVIEWERASCII);CHKERRQ(ierr);
1028*de482cd7SCarola Kruse       ierr = PetscObjectGetTabLevel((PetscObject)ilink->ksp,&tablevel);CHKERRQ(ierr);
1029*de482cd7SCarola Kruse       ierr = PetscViewerASCIISetTab(jac->gkbviewer,tablevel+1);CHKERRQ(ierr);
1030*de482cd7SCarola Kruse       ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)ilink->ksp,1);CHKERRQ(ierr);
1031*de482cd7SCarola Kruse     }
10320971522cSBarry Smith   } else {
103368bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
10340971522cSBarry Smith     i     = 0;
10350971522cSBarry Smith     ilink = jac->head;
10360971522cSBarry Smith     while (ilink) {
103723ee1639SBarry Smith       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);CHKERRQ(ierr);
10380971522cSBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
10390971522cSBarry Smith       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
10400971522cSBarry Smith       i++;
10410971522cSBarry Smith       ilink = ilink->next;
10420971522cSBarry Smith     }
10433b224e63SBarry Smith   }
10443b224e63SBarry Smith 
1045c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
10460971522cSBarry Smith   PetscFunctionReturn(0);
10470971522cSBarry Smith }
10480971522cSBarry Smith 
10495a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
1050ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
1051ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
10529df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
10535a9f2f41SSatish Balay    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
10549df09d43SBarry Smith    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
1055ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
1056ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
105779416396SBarry Smith 
10583b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
10593b224e63SBarry Smith {
10603b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10613b224e63SBarry Smith   PetscErrorCode    ierr;
10623b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
1063443836d0SMatthew G Knepley   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
10643b224e63SBarry Smith 
10653b224e63SBarry Smith   PetscFunctionBegin;
1066c5d2311dSJed Brown   switch (jac->schurfactorization) {
1067c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
1068a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
1069c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1070c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1071c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10729df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1073443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
10749df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1075c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1076c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10779df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
1078c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
10799df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
1080c096484dSStefano Zampini     ierr = VecScale(ilinkD->y,jac->schurscale);CHKERRQ(ierr);
1081c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1082c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1083c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1084c5d2311dSJed Brown     break;
1085c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
1086a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
1087c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1088c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10899df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1090443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
10919df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1092c5d2311dSJed Brown     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
1093c5d2311dSJed Brown     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
1094c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1095c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1096c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10979df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
1098c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
10999df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
1100c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1101c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1102c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1103c5d2311dSJed Brown     break;
1104c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
1105a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
1106c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1107c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11089df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
1109c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
11109df09d43SBarry 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);
1111c5d2311dSJed Brown     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
1112c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1113c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1114c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11159df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1116443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
11179df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1118c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1119c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1120c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1121c5d2311dSJed Brown     break;
1122c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
1123c238f8cdSStefano Zampini     /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1] */
11243b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11253b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11269df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1127443836d0SMatthew G Knepley     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
11289df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
11293b224e63SBarry Smith     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
11303b224e63SBarry Smith     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
11313b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11323b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11333b224e63SBarry Smith 
11349df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
11353b224e63SBarry Smith     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
11369df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
11373b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11383b224e63SBarry Smith 
1139443836d0SMatthew G Knepley     if (kspUpper == kspA) {
11403b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
11413b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
11429df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1143443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
11449df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1145443836d0SMatthew G Knepley     } else {
11469df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1147443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
11489df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1149443836d0SMatthew G Knepley       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
11509df09d43SBarry Smith       ierr = PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
1151443836d0SMatthew G Knepley       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
11529df09d43SBarry Smith       ierr = PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
1153443836d0SMatthew G Knepley       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
1154443836d0SMatthew G Knepley     }
1155c238f8cdSStefano Zampini     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11563b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11573b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1158c5d2311dSJed Brown   }
11593b224e63SBarry Smith   PetscFunctionReturn(0);
11603b224e63SBarry Smith }
11613b224e63SBarry Smith 
11620971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
11630971522cSBarry Smith {
11640971522cSBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
11650971522cSBarry Smith   PetscErrorCode     ierr;
11665a9f2f41SSatish Balay   PC_FieldSplitLink  ilink = jac->head;
1167939b8a20SBarry Smith   PetscInt           cnt,bs;
1168568ed31eSHong Zhang   KSPConvergedReason reason;
11690971522cSBarry Smith 
11700971522cSBarry Smith   PetscFunctionBegin;
117179416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
11721b9fc7fcSBarry Smith     if (jac->defaultsplit) {
1173939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
1174ce94432eSBarry 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);
1175939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
1176ce94432eSBarry 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);
11770971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
11785a9f2f41SSatish Balay       while (ilink) {
11799df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
11805a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
11819df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1182568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1183568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1184568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
1185568ed31eSHong Zhang         }
11865a9f2f41SSatish Balay         ilink = ilink->next;
11870971522cSBarry Smith       }
11880971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
11891b9fc7fcSBarry Smith     } else {
1190efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
11915a9f2f41SSatish Balay       while (ilink) {
11925a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
1193568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1194568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1195568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
1196568ed31eSHong Zhang         }
11975a9f2f41SSatish Balay         ilink = ilink->next;
11981b9fc7fcSBarry Smith       }
11991b9fc7fcSBarry Smith     }
1200e52d2c62SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
1201e52d2c62SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
1202e52d2c62SBarry Smith     /* solve on first block for first block variables */
1203e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1204e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12059df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1206e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
12079df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1208568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1209568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1210568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1211568ed31eSHong Zhang     }
1212e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1213e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1214e52d2c62SBarry Smith 
1215e52d2c62SBarry Smith     /* compute the residual only onto second block variables using first block variables */
1216e52d2c62SBarry Smith     ierr = MatMult(jac->Afield[1],ilink->y,ilink->next->x);CHKERRQ(ierr);
1217e52d2c62SBarry Smith     ilink = ilink->next;
1218e52d2c62SBarry Smith     ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
1219e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1220e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1221e52d2c62SBarry Smith 
1222e52d2c62SBarry Smith     /* solve on second block variables */
12239df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1224e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
12259df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1226568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1227568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1228568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1229568ed31eSHong Zhang     }
1230e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1231e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
123216913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
123379416396SBarry Smith     if (!jac->w1) {
123479416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
123579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
123679416396SBarry Smith     }
1237efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
12385a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
1239568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1240568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1241568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1242568ed31eSHong Zhang     }
12433e197d65SBarry Smith     cnt  = 1;
12445a9f2f41SSatish Balay     while (ilink->next) {
12455a9f2f41SSatish Balay       ilink = ilink->next;
12463e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
12473e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
12483e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
12493e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12503e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12519df09d43SBarry Smith       ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
12523e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
12539df09d43SBarry Smith       ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
12542618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
12552618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
12562618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
12572618eb8fSHong Zhang       }
12583e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12593e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12603e197d65SBarry Smith     }
126151f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
126211755939SBarry Smith       cnt -= 2;
126351f519a2SBarry Smith       while (ilink->previous) {
126451f519a2SBarry Smith         ilink = ilink->previous;
126511755939SBarry Smith         /* compute the residual only over the part of the vector needed */
126611755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
126711755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
126811755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126911755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12709df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
127111755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
12729df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1273568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1274568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1275568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
1276568ed31eSHong Zhang         }
127711755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
127811755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
127951f519a2SBarry Smith       }
128011755939SBarry Smith     }
1281ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
12820971522cSBarry Smith   PetscFunctionReturn(0);
12830971522cSBarry Smith }
12840971522cSBarry Smith 
1285a51937d4SCarola Kruse 
1286a51937d4SCarola Kruse static PetscErrorCode PCApply_FieldSplit_GKB(PC pc,Vec x,Vec y)
1287a51937d4SCarola Kruse {
1288a51937d4SCarola Kruse   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1289a51937d4SCarola Kruse   PetscErrorCode    ierr;
1290a51937d4SCarola Kruse   PC_FieldSplitLink ilinkA = jac->head,ilinkD = ilinkA->next;
1291a51937d4SCarola Kruse   KSP               ksp = ilinkA->ksp;
1292*de482cd7SCarola Kruse   Vec               u,v,Hu,d,work1,work2;
1293a51937d4SCarola Kruse   PetscScalar       alpha,z,beta,nrmz2,*vecz;
1294a51937d4SCarola Kruse   PetscReal         lowbnd,nu;
1295a51937d4SCarola Kruse   PetscInt          j,iterGKB;
1296a51937d4SCarola Kruse 
1297a51937d4SCarola Kruse   PetscFunctionBegin;
1298a51937d4SCarola Kruse   ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1299a51937d4SCarola Kruse   ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1300*de482cd7SCarola Kruse   ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1301a51937d4SCarola Kruse   ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1302a51937d4SCarola Kruse   ierr = VecDuplicate(ilinkD->x,&work1);CHKERRQ(ierr);
1303a51937d4SCarola Kruse   ierr = VecDuplicate(ilinkA->x,&work2);CHKERRQ(ierr);
1304a51937d4SCarola Kruse   ierr = VecDuplicate(ilinkD->x,&v);CHKERRQ(ierr);
1305a51937d4SCarola Kruse   ierr = VecDuplicate(ilinkD->x,&d);CHKERRQ(ierr);
1306a51937d4SCarola Kruse   ierr = VecDuplicate(ilinkA->x,&u);CHKERRQ(ierr);
1307a51937d4SCarola Kruse   ierr = VecDuplicate(ilinkA->x,&Hu);CHKERRQ(ierr);
1308a51937d4SCarola Kruse   ierr = PetscCalloc1(jac->gkbdelay,&vecz);CHKERRQ(ierr);
1309a51937d4SCarola Kruse 
1310a51937d4SCarola Kruse   /* Change RHS to comply with matrix regularization H = A + nu*B*B' */
1311a51937d4SCarola Kruse   /* Add q = q + nu*B*b */
1312a51937d4SCarola Kruse   if (jac->gkbnu) {
1313a51937d4SCarola Kruse     nu = jac->gkbnu;
1314*de482cd7SCarola Kruse     ierr = VecScale(ilinkD->x,jac->gkbnu);CHKERRQ(ierr);
1315*de482cd7SCarola Kruse     ierr = MatMultAdd(jac->B,ilinkD->x,ilinkA->x,ilinkA->x);CHKERRQ(ierr);            /* q = q + nu*B*b */
1316a51937d4SCarola Kruse   } else {
1317a51937d4SCarola Kruse     /* Situation when no augmented Lagrangian is used. Then we set inner  */
1318a51937d4SCarola Kruse     /* matrix N = I in [Ar13], and thus nu = 1.                           */
1319a51937d4SCarola Kruse     nu = 1;
1320a51937d4SCarola Kruse   }
1321a51937d4SCarola Kruse 
1322a51937d4SCarola Kruse   /* Transform rhs from [q , tilde{b}] to [0, b] */
1323*de482cd7SCarola Kruse   ierr = PetscLogEventBegin(ilinkA->event,ksp,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1324*de482cd7SCarola Kruse   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
1325*de482cd7SCarola Kruse   ierr = PetscLogEventEnd(ilinkA->event,ksp,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
1326a51937d4SCarola Kruse   ierr = MatMultHermitianTranspose(jac->B,ilinkA->y,work1);CHKERRQ(ierr);
1327*de482cd7SCarola Kruse   ierr = VecAXPBY(work1,1/nu,-1.0,ilinkD->x);CHKERRQ(ierr);             /* c = b - B'*x      */
1328a51937d4SCarola Kruse 
1329a51937d4SCarola Kruse   /* First step of algorithm */
1330*de482cd7SCarola Kruse   ierr  = VecNorm(work1,NORM_2,&beta);CHKERRQ(ierr);                  /* beta = sqrt(nu*c'*c)*/
1331*de482cd7SCarola Kruse   beta  = PetscSqrtScalar(nu)*beta;
1332*de482cd7SCarola Kruse   ierr  = VecAXPBY(v,nu/beta,0,work1);CHKERRQ(ierr);    /*SURE?*/      /* v = nu/beta *c     */
1333a51937d4SCarola Kruse   ierr  = MatMult(jac->B,v,work2);CHKERRQ(ierr);                     /* u = H^{-1}*B*v       */
1334a51937d4SCarola Kruse   ierr  = PetscLogEventBegin(ilinkA->event,ksp,work2,u,NULL);CHKERRQ(ierr);
1335a51937d4SCarola Kruse   ierr  = KSPSolve(ksp,work2,u);CHKERRQ(ierr);
1336a51937d4SCarola Kruse   ierr  = PetscLogEventEnd(ilinkA->event,ksp,work2,u,NULL);CHKERRQ(ierr);
1337a51937d4SCarola Kruse   ierr  = MatMult(jac->H,u,Hu);CHKERRQ(ierr);                      /* alpha = u'*H*u       */
1338a51937d4SCarola Kruse   ierr  = VecDot(Hu,u,&alpha);CHKERRQ(ierr);
1339*de482cd7SCarola Kruse   alpha = PetscSqrtScalar(PetscAbsScalar(alpha));
1340a51937d4SCarola Kruse   ierr  = VecScale(u,1/alpha);CHKERRQ(ierr);
1341*de482cd7SCarola Kruse   ierr  = VecAXPBY(d,1/alpha,0,v);CHKERRQ(ierr);    /*SURE?*/      /* v = nu/beta *c     */
1342*de482cd7SCarola Kruse 
1343a51937d4SCarola Kruse   z = beta/alpha;
1344a51937d4SCarola Kruse   vecz[1] = z;
1345a51937d4SCarola Kruse 
1346*de482cd7SCarola Kruse   /* Computation of first iterate x(1) and p(1) */
1347a51937d4SCarola Kruse   ierr = VecAXPY(ilinkA->y,z,u);CHKERRQ(ierr);
1348a51937d4SCarola Kruse   ierr = VecCopy(d,ilinkD->y);CHKERRQ(ierr);
1349a51937d4SCarola Kruse   ierr = VecScale(ilinkD->y,-z);CHKERRQ(ierr);
1350a51937d4SCarola Kruse 
1351a51937d4SCarola Kruse   iterGKB = 1; lowbnd = 2*jac->gkbtol;
1352*de482cd7SCarola Kruse   if (jac->gkbmonitor) {
1353*de482cd7SCarola Kruse       ierr = PetscViewerASCIIPrintf(jac->gkbviewer,"%3D GKB Lower bound estimate %14.12e\n",iterGKB,lowbnd);CHKERRQ(ierr);
1354*de482cd7SCarola Kruse   }
1355*de482cd7SCarola Kruse 
1356a51937d4SCarola Kruse   while (iterGKB < jac->gkbmaxit && lowbnd > jac->gkbtol) {
1357a51937d4SCarola Kruse     iterGKB += 1;
1358a51937d4SCarola Kruse     ierr  = MatMultHermitianTranspose(jac->B,u,work1);CHKERRQ(ierr); /* nu*(B'*u-alpha/nu*v) */
1359a51937d4SCarola Kruse     ierr  = VecScale(work1,nu);
1360a51937d4SCarola Kruse     ierr  = VecAYPX(v,-alpha,work1);CHKERRQ(ierr);
1361*de482cd7SCarola Kruse     ierr  = VecNorm(v,NORM_2,&beta);CHKERRQ(ierr);
1362*de482cd7SCarola Kruse     beta  = beta/PetscSqrtScalar(nu);
1363a51937d4SCarola Kruse     ierr  = VecScale(v,1/beta);CHKERRQ(ierr);
1364a51937d4SCarola Kruse     ierr  = MatMult(jac->B,v,work2);CHKERRQ(ierr);                 /* u = H^{-1}*B*v       */
1365a51937d4SCarola Kruse     ierr  = MatMult(jac->H,u,Hu);CHKERRQ(ierr);
1366a51937d4SCarola Kruse     ierr  = VecAXPY(work2,-beta,Hu);CHKERRQ(ierr);
1367a51937d4SCarola Kruse     ierr  = PetscLogEventBegin(ilinkA->event,ksp,work2,u,NULL);CHKERRQ(ierr);
1368a51937d4SCarola Kruse     ierr  = KSPSolve(ksp,work2,u);CHKERRQ(ierr);
1369a51937d4SCarola Kruse     ierr  = PetscLogEventEnd(ilinkA->event,ksp,work2,u,NULL);CHKERRQ(ierr);
1370a51937d4SCarola Kruse     ierr  = MatMult(jac->H,u,Hu);CHKERRQ(ierr);                    /* alpha = u'*H*u       */
1371a51937d4SCarola Kruse     ierr  = VecDot(Hu,u,&alpha);CHKERRQ(ierr);
1372*de482cd7SCarola Kruse     alpha = PetscSqrtScalar(PetscAbsScalar(alpha));
1373a51937d4SCarola Kruse     ierr  = VecScale(u,1/alpha);CHKERRQ(ierr);
1374a51937d4SCarola Kruse 
1375a51937d4SCarola Kruse     z = -beta/alpha*z;                                            /* z = beta/alpha       */
1376a51937d4SCarola Kruse     vecz[0] = z;
1377a51937d4SCarola Kruse 
1378a51937d4SCarola Kruse     /* Computation of new iterate x(i+1) and p(i+1) */
1379a51937d4SCarola Kruse     ierr = VecAXPBY(d,1/alpha,-beta/alpha,v);CHKERRQ(ierr);       /* d = 1/alpha*(v-beta*d)*/
1380a51937d4SCarola Kruse     ierr = VecAXPY(ilinkA->y,z,u);CHKERRQ(ierr);                  /* r = r + z*u           */
1381a51937d4SCarola Kruse     ierr = VecAXPY(ilinkD->y,-z,d);CHKERRQ(ierr);                 /* p = p - z*d           */
1382*de482cd7SCarola Kruse     ierr = MatMult(jac->H,ilinkA->y,Hu);CHKERRQ(ierr);            /* ||u||_H = u'*H*u      */
1383a51937d4SCarola Kruse     ierr = VecDot(Hu,ilinkA->y,&nrmz2);CHKERRQ(ierr);
1384a51937d4SCarola Kruse 
1385a51937d4SCarola Kruse     /* Compute Lower Bound estimate */
1386a51937d4SCarola Kruse     if (iterGKB > jac->gkbdelay) {
1387a51937d4SCarola Kruse       lowbnd = 0.0;
1388a51937d4SCarola Kruse       for (j=0; j<jac->gkbdelay; j++) {
1389a51937d4SCarola Kruse         lowbnd += PetscAbsScalar(vecz[j]*vecz[j]);
1390a51937d4SCarola Kruse       }
1391a51937d4SCarola Kruse       lowbnd = PetscSqrtScalar(lowbnd/PetscAbsScalar(nrmz2));
1392a51937d4SCarola Kruse     }
1393a51937d4SCarola Kruse 
1394a51937d4SCarola Kruse     for (j=0; j<jac->gkbdelay-1; j++) {
1395a51937d4SCarola Kruse       vecz[jac->gkbdelay-j-1] = vecz[jac->gkbdelay-j-2];
1396a51937d4SCarola Kruse     }
1397*de482cd7SCarola Kruse     if (jac->gkbmonitor) {
1398*de482cd7SCarola Kruse       ierr = PetscViewerASCIIPrintf(jac->gkbviewer,"%3D GKB Lower bound estimate %14.12e\n",iterGKB,lowbnd);CHKERRQ(ierr);
1399*de482cd7SCarola Kruse     }
1400a51937d4SCarola Kruse   }
1401a51937d4SCarola Kruse 
1402a51937d4SCarola Kruse   /* It would be good to have something like a gkb_monitor variable to print out the number of         */
1403a51937d4SCarola Kruse   /* iterations iterGKB and the final error estimate lowbnd. Since there is no ksp context associated  */
1404a51937d4SCarola Kruse   /* for this intermediate iteration, how can we do it?                                                */
1405a51937d4SCarola Kruse 
1406a51937d4SCarola Kruse   ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1407a51937d4SCarola Kruse   ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1408*de482cd7SCarola Kruse   ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1409a51937d4SCarola Kruse   ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1410a51937d4SCarola Kruse 
1411a51937d4SCarola Kruse   ierr = VecDestroy(&u);CHKERRQ(ierr);
1412a51937d4SCarola Kruse   ierr = VecDestroy(&v);CHKERRQ(ierr);
1413a51937d4SCarola Kruse   ierr = VecDestroy(&Hu);CHKERRQ(ierr);
1414a51937d4SCarola Kruse   ierr = VecDestroy(&d);CHKERRQ(ierr);
1415a51937d4SCarola Kruse   ierr = VecDestroy(&work1);CHKERRQ(ierr);
1416a51937d4SCarola Kruse   ierr = VecDestroy(&work2);CHKERRQ(ierr);
1417a51937d4SCarola Kruse   ierr = PetscFree(vecz);CHKERRQ(ierr);
1418a51937d4SCarola Kruse 
1419a51937d4SCarola Kruse   PetscFunctionReturn(0);
1420a51937d4SCarola Kruse }
1421a51937d4SCarola Kruse 
1422a51937d4SCarola Kruse 
1423421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1424ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1425ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
14269df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1427421e10b8SBarry Smith    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) ||                  \
14289df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1429ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1430ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
1431421e10b8SBarry Smith 
1432421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1433421e10b8SBarry Smith {
1434421e10b8SBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1435421e10b8SBarry Smith   PetscErrorCode     ierr;
1436421e10b8SBarry Smith   PC_FieldSplitLink  ilink = jac->head;
1437939b8a20SBarry Smith   PetscInt           bs;
1438291dd214SHong Zhang   KSPConvergedReason reason;
1439421e10b8SBarry Smith 
1440421e10b8SBarry Smith   PetscFunctionBegin;
1441421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1442421e10b8SBarry Smith     if (jac->defaultsplit) {
1443939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
1444ce94432eSBarry 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);
1445939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
1446ce94432eSBarry 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);
1447421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
1448421e10b8SBarry Smith       while (ilink) {
14499df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1450421e10b8SBarry Smith         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
14519df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
14522618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
14532618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
14542618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
14552618eb8fSHong Zhang         }
1456421e10b8SBarry Smith         ilink = ilink->next;
1457421e10b8SBarry Smith       }
1458421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
1459421e10b8SBarry Smith     } else {
1460421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
1461421e10b8SBarry Smith       while (ilink) {
1462421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
14632618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
14642618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
14652618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
14662618eb8fSHong Zhang         }
1467421e10b8SBarry Smith         ilink = ilink->next;
1468421e10b8SBarry Smith       }
1469421e10b8SBarry Smith     }
1470421e10b8SBarry Smith   } else {
1471421e10b8SBarry Smith     if (!jac->w1) {
1472421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
1473421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
1474421e10b8SBarry Smith     }
1475421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
1476421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1477421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
14782618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
14792618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
14802618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
14812618eb8fSHong Zhang       }
1482421e10b8SBarry Smith       while (ilink->next) {
1483421e10b8SBarry Smith         ilink = ilink->next;
14849989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1485421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1486421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1487421e10b8SBarry Smith       }
1488421e10b8SBarry Smith       while (ilink->previous) {
1489421e10b8SBarry Smith         ilink = ilink->previous;
14909989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1491421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1492421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1493421e10b8SBarry Smith       }
1494421e10b8SBarry Smith     } else {
1495421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
1496421e10b8SBarry Smith         ilink = ilink->next;
1497421e10b8SBarry Smith       }
1498421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
14992618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
15002618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
15012618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
15022618eb8fSHong Zhang       }
1503421e10b8SBarry Smith       while (ilink->previous) {
1504421e10b8SBarry Smith         ilink = ilink->previous;
15059989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1506421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1507421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1508421e10b8SBarry Smith       }
1509421e10b8SBarry Smith     }
1510421e10b8SBarry Smith   }
1511421e10b8SBarry Smith   PetscFunctionReturn(0);
1512421e10b8SBarry Smith }
1513421e10b8SBarry Smith 
1514574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
15150971522cSBarry Smith {
15160971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
15170971522cSBarry Smith   PetscErrorCode    ierr;
15185a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
15190971522cSBarry Smith 
15200971522cSBarry Smith   PetscFunctionBegin;
15215a9f2f41SSatish Balay   while (ilink) {
1522f5f0d762SBarry Smith     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1523fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
1524fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
1525443836d0SMatthew G Knepley     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
1526fcfd50ebSBarry Smith     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
1527fcfd50ebSBarry Smith     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
1528b5787286SJed Brown     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1529f5f0d762SBarry Smith     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1530f5f0d762SBarry Smith     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
1531f5f0d762SBarry Smith     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
15325a9f2f41SSatish Balay     next  = ilink->next;
1533f5f0d762SBarry Smith     ierr  = PetscFree(ilink);CHKERRQ(ierr);
15345a9f2f41SSatish Balay     ilink = next;
15350971522cSBarry Smith   }
1536f5f0d762SBarry Smith   jac->head = NULL;
153705b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1538574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
1539574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
1540574deadeSBarry Smith   } else if (jac->mat) {
15410298fd71SBarry Smith     jac->mat = NULL;
1542574deadeSBarry Smith   }
154397bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
154468dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
1545f5f0d762SBarry Smith   jac->nsplits = 0;
15466bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
15476bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
15486bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
1549a7476a74SDmitry Karpeev   ierr       = MatDestroy(&jac->schurp);CHKERRQ(ierr);
15506bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
15516bf464f9SBarry Smith   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1552d78dad28SBarry Smith   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
15536bf464f9SBarry Smith   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
15546bf464f9SBarry Smith   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
1555a51937d4SCarola Kruse   ierr       = MatDestroy(&jac->H);CHKERRQ(ierr);
1556*de482cd7SCarola Kruse   ierr       = PetscViewerDestroy(&jac->gkbviewer);CHKERRQ(ierr);
15576dbb499eSCian Wilson   jac->isrestrict = PETSC_FALSE;
1558574deadeSBarry Smith   PetscFunctionReturn(0);
1559574deadeSBarry Smith }
1560574deadeSBarry Smith 
1561574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1562574deadeSBarry Smith {
1563574deadeSBarry Smith   PetscErrorCode    ierr;
1564574deadeSBarry Smith 
1565574deadeSBarry Smith   PetscFunctionBegin;
1566574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1567c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1568285fb4e2SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurGetSubKSP_C",NULL);CHKERRQ(ierr);
1569bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr);
1570bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr);
1571bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr);
1572bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr);
1573bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr);
157429f8a81cSJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);CHKERRQ(ierr);
157537a82bf0SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);CHKERRQ(ierr);
1576bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr);
15776dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);CHKERRQ(ierr);
15780971522cSBarry Smith   PetscFunctionReturn(0);
15790971522cSBarry Smith }
15800971522cSBarry Smith 
15814416b707SBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
15820971522cSBarry Smith {
15831b9fc7fcSBarry Smith   PetscErrorCode  ierr;
15846c924f48SJed Brown   PetscInt        bs;
15857b752e3dSPatrick Sanan   PetscBool       flg;
15869dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
15873b224e63SBarry Smith   PCCompositeType ctype;
15881b9fc7fcSBarry Smith 
15890971522cSBarry Smith   PetscFunctionBegin;
1590e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"FieldSplit options");CHKERRQ(ierr);
15914ab8060aSDmitry Karpeev   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
159251f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
159351f519a2SBarry Smith   if (flg) {
159451f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
159551f519a2SBarry Smith   }
15962686e3e9SMatthew G. Knepley   jac->diag_use_amat = pc->useAmat;
15972686e3e9SMatthew 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);
15982686e3e9SMatthew G. Knepley   jac->offdiag_use_amat = pc->useAmat;
15992686e3e9SMatthew 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);
16007b752e3dSPatrick Sanan   ierr = PetscOptionsBool("-pc_fieldsplit_detect_saddle_point","Form 2-way split by detecting zero diagonal entries", "PCFieldSplitSetDetectSaddlePoint",jac->detect,&jac->detect,NULL);CHKERRQ(ierr);
16017b752e3dSPatrick Sanan   ierr = PCFieldSplitSetDetectSaddlePoint(pc,jac->detect);CHKERRQ(ierr); /* Sets split type and Schur PC type */
16023b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
16033b224e63SBarry Smith   if (flg) {
16043b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
16053b224e63SBarry Smith   }
1606c30613efSMatthew Knepley   /* Only setup fields once */
1607c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1608d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1609d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
16106c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
16116c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1612d32f9abdSBarry Smith   }
1613c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1614c5929fdfSBarry Smith     ierr = PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1615c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
16160298fd71SBarry 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);
161729f8a81cSJed 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);
1618c096484dSStefano Zampini     ierr = PetscOptionsScalar("-pc_fieldsplit_schur_scale","Scale Schur complement","PCFieldSplitSetSchurScale",jac->schurscale,&jac->schurscale,NULL);CHKERRQ(ierr);
1619a51937d4SCarola Kruse   } else if (jac->type == PC_COMPOSITE_GKB) {
1620a51937d4SCarola Kruse     ierr = PetscOptionsReal("-pc_fieldsplit_gkb_tol","The tolerance for the lower bound stopping criterion","PCFieldSplitGKBTol",jac->gkbtol,&jac->gkbtol,NULL);CHKERRQ(ierr);
1621a51937d4SCarola Kruse     ierr = PetscOptionsInt("-pc_fieldsplit_gkb_delay","The delay value for lower bound criterion","PCFieldSplitGKBDelay",jac->gkbdelay,&jac->gkbdelay,NULL);CHKERRQ(ierr);
1622a51937d4SCarola Kruse     ierr = PetscOptionsReal("-pc_fieldsplit_gkb_nu","Parameter in augmented Lagrangian approach","PCFieldSplitGKBNu",jac->gkbnu,&jac->gkbnu,NULL);CHKERRQ(ierr);
1623a51937d4SCarola Kruse     ierr = PetscOptionsInt("-pc_fieldsplit_gkb_maxit","Maximum allowed number of iterations","PCFieldSplitGKBMaxit",jac->gkbmaxit,&jac->gkbmaxit,NULL);CHKERRQ(ierr);
1624*de482cd7SCarola Kruse     ierr = PetscOptionsBool("-pc_fieldsplit_gkb_monitor","Prints number of GKB iterations and error","PCFieldSplitGKB",jac->gkbmonitor,&jac->gkbmonitor,NULL);CHKERRQ(ierr);
1625c5d2311dSJed Brown   }
16261b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
16270971522cSBarry Smith   PetscFunctionReturn(0);
16280971522cSBarry Smith }
16290971522cSBarry Smith 
16300971522cSBarry Smith /*------------------------------------------------------------------------------------*/
16310971522cSBarry Smith 
16321e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
16330971522cSBarry Smith {
163497bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
16350971522cSBarry Smith   PetscErrorCode    ierr;
16365a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
163769a612a9SBarry Smith   char              prefix[128];
16385d4c12cdSJungho Lee   PetscInt          i;
16390971522cSBarry Smith 
16400971522cSBarry Smith   PetscFunctionBegin;
16416c924f48SJed Brown   if (jac->splitdefined) {
16426c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
16436c924f48SJed Brown     PetscFunctionReturn(0);
16446c924f48SJed Brown   }
164551f519a2SBarry Smith   for (i=0; i<n; i++) {
1646e32f2f54SBarry 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);
1647e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
164851f519a2SBarry Smith   }
1649b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1650a04f6461SBarry Smith   if (splitname) {
1651db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1652a04f6461SBarry Smith   } else {
1653785e854fSJed Brown     ierr = PetscMalloc1(3,&ilink->splitname);CHKERRQ(ierr);
1654a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1655a04f6461SBarry Smith   }
16569df09d43SBarry 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 */
1657785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields);CHKERRQ(ierr);
16585a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1659785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields_col);CHKERRQ(ierr);
16605d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
16612fa5cd67SKarl Rupp 
16625a9f2f41SSatish Balay   ilink->nfields = n;
16630298fd71SBarry Smith   ilink->next    = NULL;
1664ce94432eSBarry Smith   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1665422a814eSBarry Smith   ierr           = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
166620252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
16675a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
16689005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
166969a612a9SBarry Smith 
16708caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
16715a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
16720971522cSBarry Smith 
16730971522cSBarry Smith   if (!next) {
16745a9f2f41SSatish Balay     jac->head       = ilink;
16750298fd71SBarry Smith     ilink->previous = NULL;
16760971522cSBarry Smith   } else {
16770971522cSBarry Smith     while (next->next) {
16780971522cSBarry Smith       next = next->next;
16790971522cSBarry Smith     }
16805a9f2f41SSatish Balay     next->next      = ilink;
168151f519a2SBarry Smith     ilink->previous = next;
16820971522cSBarry Smith   }
16830971522cSBarry Smith   jac->nsplits++;
16840971522cSBarry Smith   PetscFunctionReturn(0);
16850971522cSBarry Smith }
16860971522cSBarry Smith 
1687285fb4e2SStefano Zampini static PetscErrorCode  PCFieldSplitSchurGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1688285fb4e2SStefano Zampini {
1689285fb4e2SStefano Zampini   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1690285fb4e2SStefano Zampini   PetscErrorCode ierr;
1691285fb4e2SStefano Zampini 
1692285fb4e2SStefano Zampini   PetscFunctionBegin;
1693285fb4e2SStefano Zampini   *subksp = NULL;
1694285fb4e2SStefano Zampini   if (n) *n = 0;
1695285fb4e2SStefano Zampini   if (jac->type == PC_COMPOSITE_SCHUR) {
1696285fb4e2SStefano Zampini     PetscInt nn;
1697285fb4e2SStefano Zampini 
1698285fb4e2SStefano Zampini     if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitSchurGetSubKSP()");
1699285fb4e2SStefano Zampini     if (jac->nsplits != 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unexpected number of splits %D != 2",jac->nsplits);
1700285fb4e2SStefano Zampini     nn   = jac->nsplits + (jac->kspupper != jac->head->ksp ? 1 : 0);
1701285fb4e2SStefano Zampini     ierr = PetscMalloc1(nn,subksp);CHKERRQ(ierr);
1702285fb4e2SStefano Zampini     (*subksp)[0] = jac->head->ksp;
1703285fb4e2SStefano Zampini     (*subksp)[1] = jac->kspschur;
1704285fb4e2SStefano Zampini     if (jac->kspupper != jac->head->ksp) (*subksp)[2] = jac->kspupper;
1705285fb4e2SStefano Zampini     if (n) *n = nn;
1706285fb4e2SStefano Zampini   }
1707285fb4e2SStefano Zampini   PetscFunctionReturn(0);
1708285fb4e2SStefano Zampini }
1709285fb4e2SStefano Zampini 
17101e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1711e69d4d44SBarry Smith {
1712e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1713e69d4d44SBarry Smith   PetscErrorCode ierr;
1714e69d4d44SBarry Smith 
1715e69d4d44SBarry Smith   PetscFunctionBegin;
171634a3b412SBarry Smith   if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()");
1717785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
1718e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
17192fa5cd67SKarl Rupp 
1720e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
172113e0d083SBarry Smith   if (n) *n = jac->nsplits;
1722e69d4d44SBarry Smith   PetscFunctionReturn(0);
1723e69d4d44SBarry Smith }
17240971522cSBarry Smith 
17251e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
17260971522cSBarry Smith {
17270971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
17280971522cSBarry Smith   PetscErrorCode    ierr;
17290971522cSBarry Smith   PetscInt          cnt   = 0;
17305a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
17310971522cSBarry Smith 
17320971522cSBarry Smith   PetscFunctionBegin;
1733785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
17345a9f2f41SSatish Balay   while (ilink) {
17355a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
17365a9f2f41SSatish Balay     ilink            = ilink->next;
17370971522cSBarry Smith   }
17385d480477SMatthew 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);
173913e0d083SBarry Smith   if (n) *n = jac->nsplits;
17400971522cSBarry Smith   PetscFunctionReturn(0);
17410971522cSBarry Smith }
17420971522cSBarry Smith 
17436dbb499eSCian Wilson /*@C
17446dbb499eSCian Wilson     PCFieldSplitRestrictIS - Restricts the fieldsplit ISs to be within a given IS.
17456dbb499eSCian Wilson 
17466dbb499eSCian Wilson     Input Parameters:
17476dbb499eSCian Wilson +   pc  - the preconditioner context
17486dbb499eSCian Wilson +   is - the index set that defines the indices to which the fieldsplit is to be restricted
17496dbb499eSCian Wilson 
17506dbb499eSCian Wilson     Level: advanced
17516dbb499eSCian Wilson 
17526dbb499eSCian Wilson @*/
17536dbb499eSCian Wilson PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
17546dbb499eSCian Wilson {
17556dbb499eSCian Wilson   PetscErrorCode ierr;
17566dbb499eSCian Wilson 
17576dbb499eSCian Wilson   PetscFunctionBegin;
17586dbb499eSCian Wilson   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17596dbb499eSCian Wilson   PetscValidHeaderSpecific(isy,IS_CLASSID,2);
17600246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));CHKERRQ(ierr);
17616dbb499eSCian Wilson   PetscFunctionReturn(0);
17626dbb499eSCian Wilson }
17636dbb499eSCian Wilson 
17646dbb499eSCian Wilson 
17656dbb499eSCian Wilson static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
17666dbb499eSCian Wilson {
17676dbb499eSCian Wilson   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
17686dbb499eSCian Wilson   PetscErrorCode    ierr;
17696dbb499eSCian Wilson   PC_FieldSplitLink ilink = jac->head, next;
17706dbb499eSCian Wilson   PetscInt          localsize,size,sizez,i;
17716dbb499eSCian Wilson   const PetscInt    *ind, *indz;
17726dbb499eSCian Wilson   PetscInt          *indc, *indcz;
17736dbb499eSCian Wilson   PetscBool         flg;
17746dbb499eSCian Wilson 
17756dbb499eSCian Wilson   PetscFunctionBegin;
17766dbb499eSCian Wilson   ierr = ISGetLocalSize(isy,&localsize);CHKERRQ(ierr);
17776dbb499eSCian Wilson   ierr = MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));CHKERRQ(ierr);
17786dbb499eSCian Wilson   size -= localsize;
17796dbb499eSCian Wilson   while(ilink) {
17806dbb499eSCian Wilson     IS isrl,isr;
17811c7cfba8SBarry Smith     PC subpc;
17826dbb499eSCian Wilson     ierr          = ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);CHKERRQ(ierr);
17836dbb499eSCian Wilson     ierr          = ISGetLocalSize(isrl,&localsize);CHKERRQ(ierr);
17846dbb499eSCian Wilson     ierr          = PetscMalloc1(localsize,&indc);CHKERRQ(ierr);
17856dbb499eSCian Wilson     ierr          = ISGetIndices(isrl,&ind);CHKERRQ(ierr);
17866dbb499eSCian Wilson     ierr          = PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));CHKERRQ(ierr);
17876dbb499eSCian Wilson     ierr          = ISRestoreIndices(isrl,&ind);CHKERRQ(ierr);
17886dbb499eSCian Wilson     ierr          = ISDestroy(&isrl);CHKERRQ(ierr);
17896dbb499eSCian Wilson     for (i=0; i<localsize; i++) *(indc+i) += size;
17906dbb499eSCian Wilson     ierr          = ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);CHKERRQ(ierr);
17916dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
17926dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
17936dbb499eSCian Wilson     ilink->is     = isr;
17946dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
17956dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
17966dbb499eSCian Wilson     ilink->is_col = isr;
17976dbb499eSCian Wilson     ierr          = ISDestroy(&isr);CHKERRQ(ierr);
17986dbb499eSCian Wilson     ierr          = KSPGetPC(ilink->ksp, &subpc);CHKERRQ(ierr);
17996dbb499eSCian Wilson     ierr          = PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
18006dbb499eSCian Wilson     if(flg) {
18016dbb499eSCian Wilson       IS iszl,isz;
18026dbb499eSCian Wilson       MPI_Comm comm;
18036dbb499eSCian Wilson       ierr   = ISGetLocalSize(ilink->is,&localsize);CHKERRQ(ierr);
18046dbb499eSCian Wilson       comm   = PetscObjectComm((PetscObject)ilink->is);
18056dbb499eSCian Wilson       ierr   = ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);CHKERRQ(ierr);
18066dbb499eSCian Wilson       ierr   = MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
18076dbb499eSCian Wilson       sizez -= localsize;
18086dbb499eSCian Wilson       ierr   = ISGetLocalSize(iszl,&localsize);CHKERRQ(ierr);
18096dbb499eSCian Wilson       ierr   = PetscMalloc1(localsize,&indcz);CHKERRQ(ierr);
18106dbb499eSCian Wilson       ierr   = ISGetIndices(iszl,&indz);CHKERRQ(ierr);
18116dbb499eSCian Wilson       ierr   = PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));CHKERRQ(ierr);
18126dbb499eSCian Wilson       ierr   = ISRestoreIndices(iszl,&indz);CHKERRQ(ierr);
18136dbb499eSCian Wilson       ierr   = ISDestroy(&iszl);CHKERRQ(ierr);
18146dbb499eSCian Wilson       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
18156dbb499eSCian Wilson       ierr   = ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);CHKERRQ(ierr);
18166dbb499eSCian Wilson       ierr   = PCFieldSplitRestrictIS(subpc,isz);CHKERRQ(ierr);
18176dbb499eSCian Wilson       ierr   = ISDestroy(&isz);CHKERRQ(ierr);
18186dbb499eSCian Wilson     }
18196dbb499eSCian Wilson     next = ilink->next;
18206dbb499eSCian Wilson     ilink = next;
18216dbb499eSCian Wilson   }
18226dbb499eSCian Wilson   jac->isrestrict = PETSC_TRUE;
18236dbb499eSCian Wilson   PetscFunctionReturn(0);
18246dbb499eSCian Wilson }
18256dbb499eSCian Wilson 
18261e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1827704ba839SBarry Smith {
1828704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1829704ba839SBarry Smith   PetscErrorCode    ierr;
1830704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1831704ba839SBarry Smith   char              prefix[128];
1832704ba839SBarry Smith 
1833704ba839SBarry Smith   PetscFunctionBegin;
18346c924f48SJed Brown   if (jac->splitdefined) {
18356c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
18366c924f48SJed Brown     PetscFunctionReturn(0);
18376c924f48SJed Brown   }
1838b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1839a04f6461SBarry Smith   if (splitname) {
1840db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1841a04f6461SBarry Smith   } else {
1842785e854fSJed Brown     ierr = PetscMalloc1(8,&ilink->splitname);CHKERRQ(ierr);
1843b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1844a04f6461SBarry Smith   }
18459df09d43SBarry 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 */
18461ab39975SBarry Smith   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1847b5787286SJed Brown   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1848b5787286SJed Brown   ilink->is     = is;
1849b5787286SJed Brown   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1850b5787286SJed Brown   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1851b5787286SJed Brown   ilink->is_col = is;
18520298fd71SBarry Smith   ilink->next   = NULL;
1853ce94432eSBarry Smith   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1854422a814eSBarry Smith   ierr          = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
185520252d06SBarry Smith   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1856704ba839SBarry Smith   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
18579005cf84SBarry Smith   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1858704ba839SBarry Smith 
18598caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1860704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1861704ba839SBarry Smith 
1862704ba839SBarry Smith   if (!next) {
1863704ba839SBarry Smith     jac->head       = ilink;
18640298fd71SBarry Smith     ilink->previous = NULL;
1865704ba839SBarry Smith   } else {
1866704ba839SBarry Smith     while (next->next) {
1867704ba839SBarry Smith       next = next->next;
1868704ba839SBarry Smith     }
1869704ba839SBarry Smith     next->next      = ilink;
1870704ba839SBarry Smith     ilink->previous = next;
1871704ba839SBarry Smith   }
1872704ba839SBarry Smith   jac->nsplits++;
1873704ba839SBarry Smith   PetscFunctionReturn(0);
1874704ba839SBarry Smith }
1875704ba839SBarry Smith 
187694ef8ddeSSatish Balay /*@C
18770971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
18780971522cSBarry Smith 
1879ad4df100SBarry Smith     Logically Collective on PC
18800971522cSBarry Smith 
18810971522cSBarry Smith     Input Parameters:
18820971522cSBarry Smith +   pc  - the preconditioner context
18830298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
18840971522cSBarry Smith .   n - the number of fields in this split
1885db4c96c1SJed Brown -   fields - the fields in this split
18860971522cSBarry Smith 
18870971522cSBarry Smith     Level: intermediate
18880971522cSBarry Smith 
188995452b02SPatrick Sanan     Notes:
189095452b02SPatrick Sanan     Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1891d32f9abdSBarry Smith 
18927287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1893d32f9abdSBarry 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
1894d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1895d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1896d32f9abdSBarry Smith 
1897db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1898db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1899db4c96c1SJed Brown 
19005d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
19015d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
19025d4c12cdSJungho Lee      available when this routine is called.
19035d4c12cdSJungho Lee 
1904d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
19050971522cSBarry Smith 
19060971522cSBarry Smith @*/
19075d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
19080971522cSBarry Smith {
19094ac538c5SBarry Smith   PetscErrorCode ierr;
19100971522cSBarry Smith 
19110971522cSBarry Smith   PetscFunctionBegin;
19120700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1913db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1914ce94432eSBarry Smith   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1915db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
19160246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
19170971522cSBarry Smith   PetscFunctionReturn(0);
19180971522cSBarry Smith }
19190971522cSBarry Smith 
1920c84da90fSDmitry Karpeev /*@
1921c84da90fSDmitry Karpeev     PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1922c84da90fSDmitry Karpeev 
1923c84da90fSDmitry Karpeev     Logically Collective on PC
1924c84da90fSDmitry Karpeev 
1925c84da90fSDmitry Karpeev     Input Parameters:
1926c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1927c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1928c84da90fSDmitry Karpeev 
19299506b623SDmitry Karpeev     Options Database:
19309506b623SDmitry Karpeev .     -pc_fieldsplit_diag_use_amat
1931c84da90fSDmitry Karpeev 
1932c84da90fSDmitry Karpeev     Level: intermediate
1933c84da90fSDmitry Karpeev 
1934c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetDiagUseAmat(), PCFieldSplitSetOffDiagUseAmat(), PCFIELDSPLIT
1935c84da90fSDmitry Karpeev 
1936c84da90fSDmitry Karpeev @*/
1937c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1938c84da90fSDmitry Karpeev {
1939c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1940c84da90fSDmitry Karpeev   PetscBool      isfs;
1941c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1942c84da90fSDmitry Karpeev 
1943c84da90fSDmitry Karpeev   PetscFunctionBegin;
1944c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1945c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1946c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1947c84da90fSDmitry Karpeev   jac->diag_use_amat = flg;
1948c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1949c84da90fSDmitry Karpeev }
1950c84da90fSDmitry Karpeev 
1951c84da90fSDmitry Karpeev /*@
1952c84da90fSDmitry Karpeev     PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1953c84da90fSDmitry Karpeev 
1954c84da90fSDmitry Karpeev     Logically Collective on PC
1955c84da90fSDmitry Karpeev 
1956c84da90fSDmitry Karpeev     Input Parameters:
1957c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1958c84da90fSDmitry Karpeev 
1959c84da90fSDmitry Karpeev     Output Parameters:
1960c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1961c84da90fSDmitry Karpeev 
1962c84da90fSDmitry Karpeev 
1963c84da90fSDmitry Karpeev     Level: intermediate
1964c84da90fSDmitry Karpeev 
1965c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetDiagUseAmat(), PCFieldSplitGetOffDiagUseAmat(), PCFIELDSPLIT
1966c84da90fSDmitry Karpeev 
1967c84da90fSDmitry Karpeev @*/
1968c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1969c84da90fSDmitry Karpeev {
1970c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1971c84da90fSDmitry Karpeev   PetscBool      isfs;
1972c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1973c84da90fSDmitry Karpeev 
1974c84da90fSDmitry Karpeev   PetscFunctionBegin;
1975c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1976c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1977c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1978c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1979c84da90fSDmitry Karpeev   *flg = jac->diag_use_amat;
1980c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1981c84da90fSDmitry Karpeev }
1982c84da90fSDmitry Karpeev 
1983c84da90fSDmitry Karpeev /*@
1984c84da90fSDmitry Karpeev     PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1985c84da90fSDmitry Karpeev 
1986c84da90fSDmitry Karpeev     Logically Collective on PC
1987c84da90fSDmitry Karpeev 
1988c84da90fSDmitry Karpeev     Input Parameters:
1989c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1990c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1991c84da90fSDmitry Karpeev 
19929506b623SDmitry Karpeev     Options Database:
19939506b623SDmitry Karpeev .     -pc_fieldsplit_off_diag_use_amat
1994c84da90fSDmitry Karpeev 
1995c84da90fSDmitry Karpeev     Level: intermediate
1996c84da90fSDmitry Karpeev 
1997c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetOffDiagUseAmat(), PCFieldSplitSetDiagUseAmat(), PCFIELDSPLIT
1998c84da90fSDmitry Karpeev 
1999c84da90fSDmitry Karpeev @*/
2000c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
2001c84da90fSDmitry Karpeev {
2002c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2003c84da90fSDmitry Karpeev   PetscBool      isfs;
2004c84da90fSDmitry Karpeev   PetscErrorCode ierr;
2005c84da90fSDmitry Karpeev 
2006c84da90fSDmitry Karpeev   PetscFunctionBegin;
2007c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2008c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
2009c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
2010c84da90fSDmitry Karpeev   jac->offdiag_use_amat = flg;
2011c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
2012c84da90fSDmitry Karpeev }
2013c84da90fSDmitry Karpeev 
2014c84da90fSDmitry Karpeev /*@
2015c84da90fSDmitry Karpeev     PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
2016c84da90fSDmitry Karpeev 
2017c84da90fSDmitry Karpeev     Logically Collective on PC
2018c84da90fSDmitry Karpeev 
2019c84da90fSDmitry Karpeev     Input Parameters:
2020c84da90fSDmitry Karpeev .   pc  - the preconditioner object
2021c84da90fSDmitry Karpeev 
2022c84da90fSDmitry Karpeev     Output Parameters:
2023c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
2024c84da90fSDmitry Karpeev 
2025c84da90fSDmitry Karpeev 
2026c84da90fSDmitry Karpeev     Level: intermediate
2027c84da90fSDmitry Karpeev 
2028c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetOffDiagUseAmat(), PCFieldSplitGetDiagUseAmat(), PCFIELDSPLIT
2029c84da90fSDmitry Karpeev 
2030c84da90fSDmitry Karpeev @*/
2031c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
2032c84da90fSDmitry Karpeev {
2033c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2034c84da90fSDmitry Karpeev   PetscBool      isfs;
2035c84da90fSDmitry Karpeev   PetscErrorCode ierr;
2036c84da90fSDmitry Karpeev 
2037c84da90fSDmitry Karpeev   PetscFunctionBegin;
2038c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2039c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
2040c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
2041c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
2042c84da90fSDmitry Karpeev   *flg = jac->offdiag_use_amat;
2043c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
2044c84da90fSDmitry Karpeev }
2045c84da90fSDmitry Karpeev 
2046c84da90fSDmitry Karpeev 
2047c84da90fSDmitry Karpeev 
2048218d4943SBarry Smith /*@C
2049704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
2050704ba839SBarry Smith 
2051ad4df100SBarry Smith     Logically Collective on PC
2052704ba839SBarry Smith 
2053704ba839SBarry Smith     Input Parameters:
2054704ba839SBarry Smith +   pc  - the preconditioner context
20550298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
2056db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
2057704ba839SBarry Smith 
2058d32f9abdSBarry Smith 
2059a6ffb8dbSJed Brown     Notes:
2060a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
2061a6ffb8dbSJed Brown 
2062db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
2063db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
2064d32f9abdSBarry Smith 
2065704ba839SBarry Smith     Level: intermediate
2066704ba839SBarry Smith 
2067704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
2068704ba839SBarry Smith 
2069704ba839SBarry Smith @*/
20707087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
2071704ba839SBarry Smith {
20724ac538c5SBarry Smith   PetscErrorCode ierr;
2073704ba839SBarry Smith 
2074704ba839SBarry Smith   PetscFunctionBegin;
20750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
20767b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
2077db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
2078ccc738f7SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
2079704ba839SBarry Smith   PetscFunctionReturn(0);
2080704ba839SBarry Smith }
2081704ba839SBarry Smith 
208294ef8ddeSSatish Balay /*@C
208357a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
208457a9adfeSMatthew G Knepley 
208557a9adfeSMatthew G Knepley     Logically Collective on PC
208657a9adfeSMatthew G Knepley 
208757a9adfeSMatthew G Knepley     Input Parameters:
208857a9adfeSMatthew G Knepley +   pc  - the preconditioner context
208957a9adfeSMatthew G Knepley -   splitname - name of this split
209057a9adfeSMatthew G Knepley 
209157a9adfeSMatthew G Knepley     Output Parameter:
20920298fd71SBarry Smith -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
209357a9adfeSMatthew G Knepley 
209457a9adfeSMatthew G Knepley     Level: intermediate
209557a9adfeSMatthew G Knepley 
209657a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
209757a9adfeSMatthew G Knepley 
209857a9adfeSMatthew G Knepley @*/
209957a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
210057a9adfeSMatthew G Knepley {
210157a9adfeSMatthew G Knepley   PetscErrorCode ierr;
210257a9adfeSMatthew G Knepley 
210357a9adfeSMatthew G Knepley   PetscFunctionBegin;
210457a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
210557a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
210657a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
210757a9adfeSMatthew G Knepley   {
210857a9adfeSMatthew G Knepley     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
210957a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
211057a9adfeSMatthew G Knepley     PetscBool         found;
211157a9adfeSMatthew G Knepley 
21120298fd71SBarry Smith     *is = NULL;
211357a9adfeSMatthew G Knepley     while (ilink) {
211457a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
211557a9adfeSMatthew G Knepley       if (found) {
211657a9adfeSMatthew G Knepley         *is = ilink->is;
211757a9adfeSMatthew G Knepley         break;
211857a9adfeSMatthew G Knepley       }
211957a9adfeSMatthew G Knepley       ilink = ilink->next;
212057a9adfeSMatthew G Knepley     }
212157a9adfeSMatthew G Knepley   }
212257a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
212357a9adfeSMatthew G Knepley }
212457a9adfeSMatthew G Knepley 
212551f519a2SBarry Smith /*@
212651f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
212751f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
212851f519a2SBarry Smith 
2129ad4df100SBarry Smith     Logically Collective on PC
213051f519a2SBarry Smith 
213151f519a2SBarry Smith     Input Parameters:
213251f519a2SBarry Smith +   pc  - the preconditioner context
213351f519a2SBarry Smith -   bs - the block size
213451f519a2SBarry Smith 
213551f519a2SBarry Smith     Level: intermediate
213651f519a2SBarry Smith 
213751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
213851f519a2SBarry Smith 
213951f519a2SBarry Smith @*/
21407087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
214151f519a2SBarry Smith {
21424ac538c5SBarry Smith   PetscErrorCode ierr;
214351f519a2SBarry Smith 
214451f519a2SBarry Smith   PetscFunctionBegin;
21450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2146c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
21474ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
214851f519a2SBarry Smith   PetscFunctionReturn(0);
214951f519a2SBarry Smith }
215051f519a2SBarry Smith 
21510971522cSBarry Smith /*@C
215269a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
21530971522cSBarry Smith 
215469a612a9SBarry Smith    Collective on KSP
21550971522cSBarry Smith 
21560971522cSBarry Smith    Input Parameter:
21570971522cSBarry Smith .  pc - the preconditioner context
21580971522cSBarry Smith 
21590971522cSBarry Smith    Output Parameters:
216013e0d083SBarry Smith +  n - the number of splits
21613111de3cSBarry Smith -  subksp - the array of KSP contexts
21620971522cSBarry Smith 
21630971522cSBarry Smith    Note:
21649a4f7e47SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
2165d32f9abdSBarry Smith    (not the KSP just the array that contains them).
21660971522cSBarry Smith 
2167285fb4e2SStefano Zampini    You must call PCSetUp() before calling PCFieldSplitGetSubKSP().
2168285fb4e2SStefano Zampini 
2169285fb4e2SStefano Zampini    If the fieldsplit is of type PC_COMPOSITE_SCHUR, it returns the KSP object used inside the
2170285fb4e2SStefano Zampini    Schur complement and the KSP object used to iterate over the Schur complement.
2171a51937d4SCarola Kruse    To access all the KSP objects used in PC_COMPOSITE_SCHUR, use PCFieldSplitSchurGetSubKSP().
2172a51937d4SCarola Kruse 
2173a51937d4SCarola Kruse    If the fieldsplit is of type PC_COMPOSITE_GKB, it returns the KSP object used to solve the
2174a51937d4SCarola Kruse    inner linear system defined by the matrix H in each loop.
21750971522cSBarry Smith 
2176196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
21772bf68e3eSBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
2178196cc216SBarry Smith       KSP array must be.
2179196cc216SBarry Smith 
2180196cc216SBarry Smith 
21810971522cSBarry Smith    Level: advanced
21820971522cSBarry Smith 
21830971522cSBarry Smith .seealso: PCFIELDSPLIT
21840971522cSBarry Smith @*/
21857087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
21860971522cSBarry Smith {
21874ac538c5SBarry Smith   PetscErrorCode ierr;
21880971522cSBarry Smith 
21890971522cSBarry Smith   PetscFunctionBegin;
21900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
219113e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
21924ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
21930971522cSBarry Smith   PetscFunctionReturn(0);
21940971522cSBarry Smith }
21950971522cSBarry Smith 
2196285fb4e2SStefano Zampini /*@C
2197285fb4e2SStefano Zampini    PCFieldSplitSchurGetSubKSP - Gets the KSP contexts used inside the Schur complement based PCFIELDSPLIT
2198285fb4e2SStefano Zampini 
2199285fb4e2SStefano Zampini    Collective on KSP
2200285fb4e2SStefano Zampini 
2201285fb4e2SStefano Zampini    Input Parameter:
2202285fb4e2SStefano Zampini .  pc - the preconditioner context
2203285fb4e2SStefano Zampini 
2204285fb4e2SStefano Zampini    Output Parameters:
2205285fb4e2SStefano Zampini +  n - the number of splits
2206285fb4e2SStefano Zampini -  subksp - the array of KSP contexts
2207285fb4e2SStefano Zampini 
2208285fb4e2SStefano Zampini    Note:
2209285fb4e2SStefano Zampini    After PCFieldSplitSchurGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
2210285fb4e2SStefano Zampini    (not the KSP just the array that contains them).
2211285fb4e2SStefano Zampini 
2212285fb4e2SStefano Zampini    You must call PCSetUp() before calling PCFieldSplitSchurGetSubKSP().
2213285fb4e2SStefano Zampini 
2214285fb4e2SStefano Zampini    If the fieldsplit type is of type PC_COMPOSITE_SCHUR, it returns (in order)
2215285fb4e2SStefano Zampini    - the KSP used for the (1,1) block
2216285fb4e2SStefano Zampini    - the KSP used for the Schur complement (not the one used for the interior Schur solver)
2217285fb4e2SStefano Zampini    - the KSP used for the (1,1) block in the upper triangular factor (if different from that of the (1,1) block).
2218285fb4e2SStefano Zampini 
2219285fb4e2SStefano Zampini    It returns a null array if the fieldsplit is not of type PC_COMPOSITE_SCHUR; in this case, you should use PCFieldSplitGetSubKSP().
2220285fb4e2SStefano Zampini 
2221285fb4e2SStefano Zampini    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
2222285fb4e2SStefano Zampini       You can call PCFieldSplitSchurGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
2223285fb4e2SStefano Zampini       KSP array must be.
2224285fb4e2SStefano Zampini 
2225285fb4e2SStefano Zampini    Level: advanced
2226285fb4e2SStefano Zampini 
2227285fb4e2SStefano Zampini .seealso: PCFIELDSPLIT
2228285fb4e2SStefano Zampini @*/
2229285fb4e2SStefano Zampini PetscErrorCode  PCFieldSplitSchurGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
2230285fb4e2SStefano Zampini {
2231285fb4e2SStefano Zampini   PetscErrorCode ierr;
2232285fb4e2SStefano Zampini 
2233285fb4e2SStefano Zampini   PetscFunctionBegin;
2234285fb4e2SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2235285fb4e2SStefano Zampini   if (n) PetscValidIntPointer(n,2);
2236285fb4e2SStefano Zampini   ierr = PetscUseMethod(pc,"PCFieldSplitSchurGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
2237285fb4e2SStefano Zampini   PetscFunctionReturn(0);
2238285fb4e2SStefano Zampini }
2239285fb4e2SStefano Zampini 
2240e69d4d44SBarry Smith /*@
224153f2277eSBarry Smith     PCFieldSplitSetSchurPre -  Indicates what operator is used to construct the preconditioner for the Schur complement.
2242a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
2243e69d4d44SBarry Smith 
2244e69d4d44SBarry Smith     Collective on PC
2245e69d4d44SBarry Smith 
2246e69d4d44SBarry Smith     Input Parameters:
2247e69d4d44SBarry Smith +   pc      - the preconditioner context
22486fdd48a9SBarry 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
22496fdd48a9SBarry Smith               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
22500298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
2251084e4875SJed Brown 
2252e69d4d44SBarry Smith     Options Database:
22536fdd48a9SBarry Smith .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments
2254e69d4d44SBarry Smith 
2255fd1303e9SJungho Lee     Notes:
2256fd1303e9SJungho Lee $    If ptype is
2257a6a584a2SBarry Smith $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
225853f2277eSBarry Smith $             matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
2259a6a584a2SBarry Smith $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
2260a6a584a2SBarry Smith $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
2261fd1303e9SJungho Lee $             preconditioner
22626fdd48a9SBarry Smith $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
22636fdd48a9SBarry Smith $             to this function).
2264a6a584a2SBarry Smith $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
22651d71051fSDmitry Karpeev $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
22666fdd48a9SBarry Smith $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
226794350679SMatthew Knepley $        full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PCFIELDSPLIT (this is expensive)
22686fdd48a9SBarry Smith $             useful mostly as a test that the Schur complement approach can work for your problem
2269fd1303e9SJungho Lee 
2270e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
2271fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
2272a7476a74SDmitry Karpeev     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
2273fd1303e9SJungho Lee 
2274e69d4d44SBarry Smith     Level: intermediate
2275e69d4d44SBarry Smith 
22761d71051fSDmitry Karpeev .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
22771d71051fSDmitry Karpeev           MatSchurComplementSetAinvType(), PCLSC
2278e69d4d44SBarry Smith 
2279e69d4d44SBarry Smith @*/
228029f8a81cSJed Brown PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2281e69d4d44SBarry Smith {
22824ac538c5SBarry Smith   PetscErrorCode ierr;
2283e69d4d44SBarry Smith 
2284e69d4d44SBarry Smith   PetscFunctionBegin;
22850700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
228629f8a81cSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
2287e69d4d44SBarry Smith   PetscFunctionReturn(0);
2288e69d4d44SBarry Smith }
2289686bed4dSStefano Zampini 
229029f8a81cSJed Brown PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */
2291e69d4d44SBarry Smith 
229237a82bf0SJed Brown /*@
229337a82bf0SJed Brown     PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
229429f8a81cSJed Brown     preconditioned.  See PCFieldSplitSetSchurPre() for details.
229537a82bf0SJed Brown 
229637a82bf0SJed Brown     Logically Collective on PC
229737a82bf0SJed Brown 
229837a82bf0SJed Brown     Input Parameters:
229937a82bf0SJed Brown .   pc      - the preconditioner context
230037a82bf0SJed Brown 
230137a82bf0SJed Brown     Output Parameters:
230237a82bf0SJed 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
230337a82bf0SJed Brown -   userpre - matrix to use for preconditioning (with PC_FIELDSPLIT_PRE_USER), or NULL
230437a82bf0SJed Brown 
230537a82bf0SJed Brown     Level: intermediate
230637a82bf0SJed Brown 
230729f8a81cSJed Brown .seealso: PCFieldSplitSetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
230837a82bf0SJed Brown 
230937a82bf0SJed Brown @*/
231037a82bf0SJed Brown PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
231137a82bf0SJed Brown {
231237a82bf0SJed Brown   PetscErrorCode ierr;
231337a82bf0SJed Brown 
231437a82bf0SJed Brown   PetscFunctionBegin;
231537a82bf0SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
231637a82bf0SJed Brown   ierr = PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));CHKERRQ(ierr);
2317e69d4d44SBarry Smith   PetscFunctionReturn(0);
2318e69d4d44SBarry Smith }
2319e69d4d44SBarry Smith 
232045e7fc46SDmitry Karpeev /*@
2321470b340bSDmitry Karpeev     PCFieldSplitSchurGetS -  extract the MatSchurComplement object used by this PC in case it needs to be configured separately
232245e7fc46SDmitry Karpeev 
232345e7fc46SDmitry Karpeev     Not collective
232445e7fc46SDmitry Karpeev 
232545e7fc46SDmitry Karpeev     Input Parameter:
232645e7fc46SDmitry Karpeev .   pc      - the preconditioner context
232745e7fc46SDmitry Karpeev 
232845e7fc46SDmitry Karpeev     Output Parameter:
2329470b340bSDmitry Karpeev .   S       - the Schur complement matrix
233045e7fc46SDmitry Karpeev 
2331470b340bSDmitry Karpeev     Notes:
2332470b340bSDmitry Karpeev     This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().
233345e7fc46SDmitry Karpeev 
233445e7fc46SDmitry Karpeev     Level: advanced
233545e7fc46SDmitry Karpeev 
233629f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurRestoreS()
233745e7fc46SDmitry Karpeev 
233845e7fc46SDmitry Karpeev @*/
2339470b340bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
234045e7fc46SDmitry Karpeev {
234145e7fc46SDmitry Karpeev   PetscErrorCode ierr;
234245e7fc46SDmitry Karpeev   const char*    t;
234345e7fc46SDmitry Karpeev   PetscBool      isfs;
234445e7fc46SDmitry Karpeev   PC_FieldSplit  *jac;
234545e7fc46SDmitry Karpeev 
234645e7fc46SDmitry Karpeev   PetscFunctionBegin;
234745e7fc46SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
234845e7fc46SDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
234945e7fc46SDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
235045e7fc46SDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
235145e7fc46SDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
2352470b340bSDmitry 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);
2353470b340bSDmitry Karpeev   if (S) *S = jac->schur;
235445e7fc46SDmitry Karpeev   PetscFunctionReturn(0);
235545e7fc46SDmitry Karpeev }
235645e7fc46SDmitry Karpeev 
2357470b340bSDmitry Karpeev /*@
2358470b340bSDmitry Karpeev     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC
2359470b340bSDmitry Karpeev 
2360470b340bSDmitry Karpeev     Not collective
2361470b340bSDmitry Karpeev 
2362470b340bSDmitry Karpeev     Input Parameters:
2363470b340bSDmitry Karpeev +   pc      - the preconditioner context
2364470b340bSDmitry Karpeev .   S       - the Schur complement matrix
2365470b340bSDmitry Karpeev 
2366470b340bSDmitry Karpeev     Level: advanced
2367470b340bSDmitry Karpeev 
236829f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurGetS()
2369470b340bSDmitry Karpeev 
2370470b340bSDmitry Karpeev @*/
2371bca69d2bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
2372470b340bSDmitry Karpeev {
2373470b340bSDmitry Karpeev   PetscErrorCode ierr;
2374470b340bSDmitry Karpeev   const char*    t;
2375470b340bSDmitry Karpeev   PetscBool      isfs;
2376470b340bSDmitry Karpeev   PC_FieldSplit  *jac;
2377470b340bSDmitry Karpeev 
2378470b340bSDmitry Karpeev   PetscFunctionBegin;
2379470b340bSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2380470b340bSDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
2381470b340bSDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
2382470b340bSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2383470b340bSDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
2384470b340bSDmitry 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);
2385bca69d2bSDmitry Karpeev   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
2386470b340bSDmitry Karpeev   PetscFunctionReturn(0);
2387470b340bSDmitry Karpeev }
2388470b340bSDmitry Karpeev 
2389470b340bSDmitry Karpeev 
239029f8a81cSJed Brown static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2391e69d4d44SBarry Smith {
2392e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2393084e4875SJed Brown   PetscErrorCode ierr;
2394e69d4d44SBarry Smith 
2395e69d4d44SBarry Smith   PetscFunctionBegin;
2396084e4875SJed Brown   jac->schurpre = ptype;
2397a7476a74SDmitry Karpeev   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
23986bf464f9SBarry Smith     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
2399084e4875SJed Brown     jac->schur_user = pre;
2400084e4875SJed Brown     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
2401084e4875SJed Brown   }
2402e69d4d44SBarry Smith   PetscFunctionReturn(0);
2403e69d4d44SBarry Smith }
2404e69d4d44SBarry Smith 
240537a82bf0SJed Brown static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
240637a82bf0SJed Brown {
240737a82bf0SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
240837a82bf0SJed Brown 
240937a82bf0SJed Brown   PetscFunctionBegin;
241037a82bf0SJed Brown   *ptype = jac->schurpre;
241137a82bf0SJed Brown   *pre   = jac->schur_user;
241237a82bf0SJed Brown   PetscFunctionReturn(0);
241337a82bf0SJed Brown }
241437a82bf0SJed Brown 
2415ab1df9f5SJed Brown /*@
24160ffb0e17SBarry Smith     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain in the preconditioner
2417ab1df9f5SJed Brown 
2418ab1df9f5SJed Brown     Collective on PC
2419ab1df9f5SJed Brown 
2420ab1df9f5SJed Brown     Input Parameters:
2421ab1df9f5SJed Brown +   pc  - the preconditioner context
2422c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
2423ab1df9f5SJed Brown 
2424ab1df9f5SJed Brown     Options Database:
2425c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
2426ab1df9f5SJed Brown 
2427ab1df9f5SJed Brown 
2428ab1df9f5SJed Brown     Level: intermediate
2429ab1df9f5SJed Brown 
2430ab1df9f5SJed Brown     Notes:
2431ab1df9f5SJed Brown     The FULL factorization is
2432ab1df9f5SJed Brown 
24330ffb0e17SBarry Smith $   (A   B)  = (1       0) (A   0) (1  Ainv*B)  = L D U
24340ffb0e17SBarry Smith $   (C   E)    (C*Ainv  1) (0   S) (0     1  )
2435ab1df9f5SJed Brown 
24360ffb0e17SBarry 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,
2437c096484dSStefano 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().
2438ab1df9f5SJed Brown 
24390ffb0e17SBarry Smith $    If A and S are solved exactly
24400ffb0e17SBarry Smith $      *) FULL factorization is a direct solver.
24410ffb0e17SBarry 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.
24420ffb0e17SBarry 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.
2443ab1df9f5SJed Brown 
24440ffb0e17SBarry Smith     If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
24450ffb0e17SBarry 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.
24460ffb0e17SBarry Smith 
24470ffb0e17SBarry Smith     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES.
24480ffb0e17SBarry Smith 
24490ffb0e17SBarry 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).
2450ab1df9f5SJed Brown 
2451ab1df9f5SJed Brown     References:
245296a0c994SBarry Smith +   1. - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
245396a0c994SBarry Smith -   2. - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).
2454ab1df9f5SJed Brown 
2455c096484dSStefano Zampini .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCFieldSplitSetSchurScale()
2456ab1df9f5SJed Brown @*/
2457c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2458ab1df9f5SJed Brown {
2459ab1df9f5SJed Brown   PetscErrorCode ierr;
2460ab1df9f5SJed Brown 
2461ab1df9f5SJed Brown   PetscFunctionBegin;
2462ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2463c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
2464ab1df9f5SJed Brown   PetscFunctionReturn(0);
2465ab1df9f5SJed Brown }
2466ab1df9f5SJed Brown 
24671e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2468ab1df9f5SJed Brown {
2469ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2470ab1df9f5SJed Brown 
2471ab1df9f5SJed Brown   PetscFunctionBegin;
2472ab1df9f5SJed Brown   jac->schurfactorization = ftype;
2473ab1df9f5SJed Brown   PetscFunctionReturn(0);
2474ab1df9f5SJed Brown }
2475ab1df9f5SJed Brown 
2476c096484dSStefano Zampini /*@
2477c096484dSStefano Zampini     PCFieldSplitSetSchurScale -  Controls the sign flip of S for PC_FIELDSPLIT_SCHUR_FACT_DIAG.
2478c096484dSStefano Zampini 
2479c096484dSStefano Zampini     Collective on PC
2480c096484dSStefano Zampini 
2481c096484dSStefano Zampini     Input Parameters:
2482c096484dSStefano Zampini +   pc    - the preconditioner context
2483c096484dSStefano Zampini -   scale - scaling factor for the Schur complement
2484c096484dSStefano Zampini 
2485c096484dSStefano Zampini     Options Database:
2486c096484dSStefano Zampini .     -pc_fieldsplit_schur_scale - default is -1.0
2487c096484dSStefano Zampini 
2488c096484dSStefano Zampini     Level: intermediate
2489c096484dSStefano Zampini 
2490c096484dSStefano Zampini .seealso: PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurFactType, PCFieldSplitSetSchurScale()
2491c096484dSStefano Zampini @*/
2492c096484dSStefano Zampini PetscErrorCode PCFieldSplitSetSchurScale(PC pc,PetscScalar scale)
2493c096484dSStefano Zampini {
2494c096484dSStefano Zampini   PetscErrorCode ierr;
2495c096484dSStefano Zampini 
2496c096484dSStefano Zampini   PetscFunctionBegin;
2497c096484dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2498c096484dSStefano Zampini   PetscValidLogicalCollectiveScalar(pc,scale,2);
2499c096484dSStefano Zampini   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurScale_C",(PC,PetscScalar),(pc,scale));CHKERRQ(ierr);
2500c096484dSStefano Zampini   PetscFunctionReturn(0);
2501c096484dSStefano Zampini }
2502c096484dSStefano Zampini 
2503c096484dSStefano Zampini static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc,PetscScalar scale)
2504c096484dSStefano Zampini {
2505c096484dSStefano Zampini   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2506c096484dSStefano Zampini 
2507c096484dSStefano Zampini   PetscFunctionBegin;
2508c096484dSStefano Zampini   jac->schurscale = scale;
2509c096484dSStefano Zampini   PetscFunctionReturn(0);
2510c096484dSStefano Zampini }
2511c096484dSStefano Zampini 
251230ad9308SMatthew Knepley /*@C
25138c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
251430ad9308SMatthew Knepley 
251530ad9308SMatthew Knepley    Collective on KSP
251630ad9308SMatthew Knepley 
251730ad9308SMatthew Knepley    Input Parameter:
251830ad9308SMatthew Knepley .  pc - the preconditioner context
251930ad9308SMatthew Knepley 
252030ad9308SMatthew Knepley    Output Parameters:
2521a04f6461SBarry Smith +  A00 - the (0,0) block
2522a04f6461SBarry Smith .  A01 - the (0,1) block
2523a04f6461SBarry Smith .  A10 - the (1,0) block
2524a04f6461SBarry Smith -  A11 - the (1,1) block
252530ad9308SMatthew Knepley 
252630ad9308SMatthew Knepley    Level: advanced
252730ad9308SMatthew Knepley 
252830ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
252930ad9308SMatthew Knepley @*/
2530a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
253130ad9308SMatthew Knepley {
253230ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
253330ad9308SMatthew Knepley 
253430ad9308SMatthew Knepley   PetscFunctionBegin;
25350700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2536ce94432eSBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2537a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
2538a04f6461SBarry Smith   if (A01) *A01 = jac->B;
2539a04f6461SBarry Smith   if (A10) *A10 = jac->C;
2540a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
254130ad9308SMatthew Knepley   PetscFunctionReturn(0);
254230ad9308SMatthew Knepley }
254330ad9308SMatthew Knepley 
2544b09e027eSCarola Kruse /*@
2545b09e027eSCarola Kruse     PCFieldSplitSetGKBTol -  Sets the solver tolerance for the Golub-Kahan bidiagonalization preconditioner.
2546b09e027eSCarola Kruse 
2547b09e027eSCarola Kruse     Collective on PC
2548b09e027eSCarola Kruse 
2549b09e027eSCarola Kruse     Input Parameters:
2550b09e027eSCarola Kruse +   pc        - the preconditioner context
2551b09e027eSCarola Kruse -   tolerance - the solver tolerance
2552b09e027eSCarola Kruse 
2553b09e027eSCarola Kruse     Options Database:
2554b09e027eSCarola Kruse .     -pc_fieldsplit_gkb_tol - default is 1e-5
2555b09e027eSCarola Kruse 
2556b09e027eSCarola Kruse     Level: intermediate
2557b09e027eSCarola Kruse 
2558b09e027eSCarola Kruse .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBNu(), PCFieldSplitSetGKBMaxit()
2559b09e027eSCarola Kruse @*/
2560b09e027eSCarola Kruse PetscErrorCode PCFieldSplitSetGKBTol(PC pc,PetscReal tolerance)
2561b09e027eSCarola Kruse {
2562b09e027eSCarola Kruse   PetscErrorCode ierr;
2563b09e027eSCarola Kruse 
2564b09e027eSCarola Kruse   PetscFunctionBegin;
2565b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2566*de482cd7SCarola Kruse   PetscValidLogicalCollectiveReal(pc,tolerance,2);
2567b09e027eSCarola Kruse   ierr = PetscTryMethod(pc,"PCFieldSplitSetGKBTol_C",(PC,PetscReal),(pc,tolerance));CHKERRQ(ierr);
2568b09e027eSCarola Kruse   PetscFunctionReturn(0);
2569b09e027eSCarola Kruse }
2570b09e027eSCarola Kruse 
2571b09e027eSCarola Kruse static PetscErrorCode PCFieldSplitSetGKBTol_FieldSplit(PC pc,PetscReal tolerance)
2572b09e027eSCarola Kruse {
2573b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2574b09e027eSCarola Kruse 
2575b09e027eSCarola Kruse   PetscFunctionBegin;
2576b09e027eSCarola Kruse   jac->gkbtol = tolerance;
2577b09e027eSCarola Kruse   PetscFunctionReturn(0);
2578b09e027eSCarola Kruse }
2579b09e027eSCarola Kruse 
2580b09e027eSCarola Kruse 
2581b09e027eSCarola Kruse /*@
2582b09e027eSCarola Kruse     PCFieldSplitSetGKBMaxit -  Sets the maximum number of iterations for the Golub-Kahan bidiagonalization preconditioner.
2583b09e027eSCarola Kruse 
2584b09e027eSCarola Kruse     Collective on PC
2585b09e027eSCarola Kruse 
2586b09e027eSCarola Kruse     Input Parameters:
2587b09e027eSCarola Kruse +   pc     - the preconditioner context
2588b09e027eSCarola Kruse -   maxit  - the maximum number of iterations
2589b09e027eSCarola Kruse 
2590b09e027eSCarola Kruse     Options Database:
2591b09e027eSCarola Kruse .     -pc_fieldsplit_gkb_maxit - default is 100
2592b09e027eSCarola Kruse 
2593b09e027eSCarola Kruse     Level: intermediate
2594b09e027eSCarola Kruse 
2595b09e027eSCarola Kruse .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBNu()
2596b09e027eSCarola Kruse @*/
2597b09e027eSCarola Kruse PetscErrorCode PCFieldSplitSetGKBMaxit(PC pc,PetscInt maxit)
2598b09e027eSCarola Kruse {
2599b09e027eSCarola Kruse   PetscErrorCode ierr;
2600b09e027eSCarola Kruse 
2601b09e027eSCarola Kruse   PetscFunctionBegin;
2602b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2603*de482cd7SCarola Kruse   PetscValidLogicalCollectiveInt(pc,maxit,2);
2604b09e027eSCarola Kruse   ierr = PetscTryMethod(pc,"PCFieldSplitSetGKBMaxit_C",(PC,PetscInt),(pc,maxit));CHKERRQ(ierr);
2605b09e027eSCarola Kruse   PetscFunctionReturn(0);
2606b09e027eSCarola Kruse }
2607b09e027eSCarola Kruse 
2608b09e027eSCarola Kruse static PetscErrorCode PCFieldSplitSetGKBMaxit_FieldSplit(PC pc,PetscInt maxit)
2609b09e027eSCarola Kruse {
2610b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2611b09e027eSCarola Kruse 
2612b09e027eSCarola Kruse   PetscFunctionBegin;
2613b09e027eSCarola Kruse   jac->gkbmaxit = maxit;
2614b09e027eSCarola Kruse   PetscFunctionReturn(0);
2615b09e027eSCarola Kruse }
2616b09e027eSCarola Kruse 
2617b09e027eSCarola Kruse /*@
2618b09e027eSCarola Kruse     PCFieldSplitSetGKBDelay -  Sets the delay in the lower bound error estimate in the Golub-Kahan bidiagonalization preconditioner.
2619b09e027eSCarola Kruse 
2620b09e027eSCarola Kruse     Collective on PC
2621b09e027eSCarola Kruse 
2622b09e027eSCarola Kruse     Notes:
2623b09e027eSCarola Kruse     The algorithm uses a lower bound estimate of the error in energy norm as stopping criterion. For more details on the
2624b09e027eSCarola Kruse     Golub-Kahan bidiagonalization method and its lower bound stopping criterion, please refer to
2625b09e027eSCarola Kruse 
2626b09e027eSCarola Kruse [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2627b09e027eSCarola Kruse 
2628b09e027eSCarola Kruse     Input Parameters:
2629b09e027eSCarola Kruse +   pc     - the preconditioner context
2630b09e027eSCarola Kruse -   delay  - the delay window in the lower bound estimate
2631b09e027eSCarola Kruse 
2632b09e027eSCarola Kruse     Options Database:
2633b09e027eSCarola Kruse .     -pc_fieldsplit_gkb_delay - default is 5
2634b09e027eSCarola Kruse 
2635b09e027eSCarola Kruse     Level: intermediate
2636b09e027eSCarola Kruse 
2637b09e027eSCarola Kruse .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBNu(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBMaxit()
2638b09e027eSCarola Kruse @*/
2639b09e027eSCarola Kruse PetscErrorCode PCFieldSplitSetGKBDelay(PC pc,PetscInt delay)
2640b09e027eSCarola Kruse {
2641b09e027eSCarola Kruse   PetscErrorCode ierr;
2642b09e027eSCarola Kruse 
2643b09e027eSCarola Kruse   PetscFunctionBegin;
2644b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2645b09e027eSCarola Kruse   PetscValidLogicalCollectiveInt(pc,delay,2);
2646b09e027eSCarola Kruse   ierr = PetscTryMethod(pc,"PCFieldSplitSetGKBDelay_C",(PC,PetscInt),(pc,delay));CHKERRQ(ierr);
2647b09e027eSCarola Kruse   PetscFunctionReturn(0);
2648b09e027eSCarola Kruse }
2649b09e027eSCarola Kruse 
2650b09e027eSCarola Kruse static PetscErrorCode PCFieldSplitSetGKBDelay_FieldSplit(PC pc,PetscInt delay)
2651b09e027eSCarola Kruse {
2652b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2653b09e027eSCarola Kruse 
2654b09e027eSCarola Kruse   PetscFunctionBegin;
2655b09e027eSCarola Kruse   jac->gkbdelay = delay;
2656b09e027eSCarola Kruse   PetscFunctionReturn(0);
2657b09e027eSCarola Kruse }
2658b09e027eSCarola Kruse 
2659b09e027eSCarola Kruse /*@
2660b09e027eSCarola Kruse     PCFieldSplitSetGKBNu -  Sets the scalar value nu >= 0 in the transformation H = A00 + nu*A01*A01' of the (1,1) block in the Golub-Kahan bidiagonalization preconditioner.
2661b09e027eSCarola Kruse 
2662b09e027eSCarola Kruse     Collective on PC
2663b09e027eSCarola Kruse 
2664b09e027eSCarola Kruse     Notes:
2665b09e027eSCarola Kruse     This shift is in general done to obtain better convergence properties for the outer loop in the algorithm. This is often achieved by chosing nu sufficiently large. However,
2666b09e027eSCarola Kruse     if nu is chosen too big, the matrix H might be badly conditioned and the solution of the linear system Hx = b in the inner loop gets difficult. It is therefore
2667b09e027eSCarola Kruse     necessary to find a good balance in between the convergence of the inner and outer loop.
2668b09e027eSCarola Kruse 
2669b09e027eSCarola Kruse     For nu = 0, no shift is done. In this case A00 has to be positive definite. The matrix N in [Ar13] is then chosen as identity.
2670b09e027eSCarola Kruse 
2671b09e027eSCarola Kruse [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2672b09e027eSCarola Kruse 
2673b09e027eSCarola Kruse     Input Parameters:
2674b09e027eSCarola Kruse +   pc     - the preconditioner context
2675b09e027eSCarola Kruse -   nu     - the shift parameter
2676b09e027eSCarola Kruse 
2677b09e027eSCarola Kruse     Options Database:
2678b09e027eSCarola Kruse .     -pc_fieldsplit_gkb_nu - default is 1
2679b09e027eSCarola Kruse 
2680b09e027eSCarola Kruse     Level: intermediate
2681b09e027eSCarola Kruse 
2682b09e027eSCarola Kruse .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBMaxit()
2683b09e027eSCarola Kruse @*/
2684b09e027eSCarola Kruse PetscErrorCode PCFieldSplitSetGKBNu(PC pc,PetscReal nu)
2685b09e027eSCarola Kruse {
2686b09e027eSCarola Kruse   PetscErrorCode ierr;
2687b09e027eSCarola Kruse 
2688b09e027eSCarola Kruse   PetscFunctionBegin;
2689b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2690*de482cd7SCarola Kruse   PetscValidLogicalCollectiveReal(pc,nu,2);
2691b09e027eSCarola Kruse   ierr = PetscTryMethod(pc,"PCFieldSplitSetGKBNu_C",(PC,PetscReal),(pc,nu));CHKERRQ(ierr);
2692b09e027eSCarola Kruse   PetscFunctionReturn(0);
2693b09e027eSCarola Kruse }
2694b09e027eSCarola Kruse 
2695b09e027eSCarola Kruse static PetscErrorCode PCFieldSplitSetGKBNu_FieldSplit(PC pc,PetscReal nu)
2696b09e027eSCarola Kruse {
2697b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2698b09e027eSCarola Kruse 
2699b09e027eSCarola Kruse   PetscFunctionBegin;
2700b09e027eSCarola Kruse   jac->gkbnu = nu;
2701b09e027eSCarola Kruse   PetscFunctionReturn(0);
2702b09e027eSCarola Kruse }
2703b09e027eSCarola Kruse 
2704b09e027eSCarola Kruse 
27051e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
270679416396SBarry Smith {
270779416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2708e69d4d44SBarry Smith   PetscErrorCode ierr;
270979416396SBarry Smith 
271079416396SBarry Smith   PetscFunctionBegin;
271179416396SBarry Smith   jac->type = type;
27123b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
27133b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
27143b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
27152fa5cd67SKarl Rupp 
2716bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
271729f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);CHKERRQ(ierr);
271837a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);CHKERRQ(ierr);
2719bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
2720c096484dSStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",PCFieldSplitSetSchurScale_FieldSplit);CHKERRQ(ierr);
2721b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBTol_C",0);CHKERRQ(ierr);
2722b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBMaxit_C",0);CHKERRQ(ierr);
2723b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBNu_C",0);CHKERRQ(ierr);
2724b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBDelay_C",0);CHKERRQ(ierr);
2725a51937d4SCarola Kruse   } else if (type == PC_COMPOSITE_GKB){
2726a51937d4SCarola Kruse     pc->ops->apply = PCApply_FieldSplit_GKB;
2727a51937d4SCarola Kruse     pc->ops->view  = PCView_FieldSplit_GKB;
2728e69d4d44SBarry Smith 
2729a51937d4SCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
2730a51937d4SCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);CHKERRQ(ierr);
2731a51937d4SCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);CHKERRQ(ierr);
2732a51937d4SCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
2733a51937d4SCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",0);CHKERRQ(ierr);
2734b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBTol_C",PCFieldSplitSetGKBTol_FieldSplit);CHKERRQ(ierr);
2735b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBMaxit_C",PCFieldSplitSetGKBMaxit_FieldSplit);CHKERRQ(ierr);
2736b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBNu_C",PCFieldSplitSetGKBNu_FieldSplit);CHKERRQ(ierr);
2737b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBDelay_C",PCFieldSplitSetGKBDelay_FieldSplit);CHKERRQ(ierr);
27383b224e63SBarry Smith   } else {
27393b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
27403b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
27412fa5cd67SKarl Rupp 
2742bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
274329f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);CHKERRQ(ierr);
274437a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);CHKERRQ(ierr);
2745bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
2746c096484dSStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",0);CHKERRQ(ierr);
2747b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBTol_C",0);CHKERRQ(ierr);
2748b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBMaxit_C",0);CHKERRQ(ierr);
2749b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBNu_C",0);CHKERRQ(ierr);
2750b09e027eSCarola Kruse     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBDelay_C",0);CHKERRQ(ierr);
27513b224e63SBarry Smith   }
275279416396SBarry Smith   PetscFunctionReturn(0);
275379416396SBarry Smith }
275479416396SBarry Smith 
27551e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
275651f519a2SBarry Smith {
275751f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
275851f519a2SBarry Smith 
275951f519a2SBarry Smith   PetscFunctionBegin;
2760ce94432eSBarry Smith   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2761ce94432eSBarry 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);
276251f519a2SBarry Smith   jac->bs = bs;
276351f519a2SBarry Smith   PetscFunctionReturn(0);
276451f519a2SBarry Smith }
276551f519a2SBarry Smith 
2766bc08b0f1SBarry Smith /*@
276779416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
276879416396SBarry Smith 
276979416396SBarry Smith    Collective on PC
277079416396SBarry Smith 
277179416396SBarry Smith    Input Parameter:
277279416396SBarry Smith .  pc - the preconditioner context
277381540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
277479416396SBarry Smith 
277579416396SBarry Smith    Options Database Key:
2776a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
277779416396SBarry Smith 
2778b02e2d75SMatthew G Knepley    Level: Intermediate
277979416396SBarry Smith 
278079416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
278179416396SBarry Smith 
278279416396SBarry Smith .seealso: PCCompositeSetType()
278379416396SBarry Smith 
278479416396SBarry Smith @*/
27857087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
278679416396SBarry Smith {
27874ac538c5SBarry Smith   PetscErrorCode ierr;
278879416396SBarry Smith 
278979416396SBarry Smith   PetscFunctionBegin;
27900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
27914ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
279279416396SBarry Smith   PetscFunctionReturn(0);
279379416396SBarry Smith }
279479416396SBarry Smith 
2795b02e2d75SMatthew G Knepley /*@
2796b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
2797b02e2d75SMatthew G Knepley 
2798b02e2d75SMatthew G Knepley   Not collective
2799b02e2d75SMatthew G Knepley 
2800b02e2d75SMatthew G Knepley   Input Parameter:
2801b02e2d75SMatthew G Knepley . pc - the preconditioner context
2802b02e2d75SMatthew G Knepley 
2803b02e2d75SMatthew G Knepley   Output Parameter:
2804b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
2805b02e2d75SMatthew G Knepley 
2806b02e2d75SMatthew G Knepley   Level: Intermediate
2807b02e2d75SMatthew G Knepley 
2808b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2809b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
2810b02e2d75SMatthew G Knepley @*/
2811b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2812b02e2d75SMatthew G Knepley {
2813b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2814b02e2d75SMatthew G Knepley 
2815b02e2d75SMatthew G Knepley   PetscFunctionBegin;
2816b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2817b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
2818b02e2d75SMatthew G Knepley   *type = jac->type;
2819b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
2820b02e2d75SMatthew G Knepley }
2821b02e2d75SMatthew G Knepley 
28224ab8060aSDmitry Karpeev /*@
28234ab8060aSDmitry Karpeev    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
28244ab8060aSDmitry Karpeev 
28254ab8060aSDmitry Karpeev    Logically Collective
28264ab8060aSDmitry Karpeev 
28274ab8060aSDmitry Karpeev    Input Parameters:
28284ab8060aSDmitry Karpeev +  pc   - the preconditioner context
28294ab8060aSDmitry Karpeev -  flg  - boolean indicating whether to use field splits defined by the DM
28304ab8060aSDmitry Karpeev 
28314ab8060aSDmitry Karpeev    Options Database Key:
28324ab8060aSDmitry Karpeev .  -pc_fieldsplit_dm_splits
28334ab8060aSDmitry Karpeev 
28344ab8060aSDmitry Karpeev    Level: Intermediate
28354ab8060aSDmitry Karpeev 
28364ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
28374ab8060aSDmitry Karpeev 
28384ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits()
28394ab8060aSDmitry Karpeev 
28404ab8060aSDmitry Karpeev @*/
28414ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
28424ab8060aSDmitry Karpeev {
28434ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
28444ab8060aSDmitry Karpeev   PetscBool      isfs;
28454ab8060aSDmitry Karpeev   PetscErrorCode ierr;
28464ab8060aSDmitry Karpeev 
28474ab8060aSDmitry Karpeev   PetscFunctionBegin;
28484ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
28494ab8060aSDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
28504ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
28514ab8060aSDmitry Karpeev   if (isfs) {
28524ab8060aSDmitry Karpeev     jac->dm_splits = flg;
28534ab8060aSDmitry Karpeev   }
28544ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
28554ab8060aSDmitry Karpeev }
28564ab8060aSDmitry Karpeev 
28574ab8060aSDmitry Karpeev 
28584ab8060aSDmitry Karpeev /*@
28594ab8060aSDmitry Karpeev    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
28604ab8060aSDmitry Karpeev 
28614ab8060aSDmitry Karpeev    Logically Collective
28624ab8060aSDmitry Karpeev 
28634ab8060aSDmitry Karpeev    Input Parameter:
28644ab8060aSDmitry Karpeev .  pc   - the preconditioner context
28654ab8060aSDmitry Karpeev 
28664ab8060aSDmitry Karpeev    Output Parameter:
28674ab8060aSDmitry Karpeev .  flg  - boolean indicating whether to use field splits defined by the DM
28684ab8060aSDmitry Karpeev 
28694ab8060aSDmitry Karpeev    Level: Intermediate
28704ab8060aSDmitry Karpeev 
28714ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
28724ab8060aSDmitry Karpeev 
28734ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits()
28744ab8060aSDmitry Karpeev 
28754ab8060aSDmitry Karpeev @*/
28764ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
28774ab8060aSDmitry Karpeev {
28784ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
28794ab8060aSDmitry Karpeev   PetscBool      isfs;
28804ab8060aSDmitry Karpeev   PetscErrorCode ierr;
28814ab8060aSDmitry Karpeev 
28824ab8060aSDmitry Karpeev   PetscFunctionBegin;
28834ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
28844ab8060aSDmitry Karpeev   PetscValidPointer(flg,2);
28854ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
28864ab8060aSDmitry Karpeev   if (isfs) {
28874ab8060aSDmitry Karpeev     if(flg) *flg = jac->dm_splits;
28884ab8060aSDmitry Karpeev   }
28894ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
28904ab8060aSDmitry Karpeev }
28914ab8060aSDmitry Karpeev 
28927b752e3dSPatrick Sanan /*@
28937b752e3dSPatrick Sanan    PCFieldSplitGetDetectSaddlePoint - Returns flag indicating whether PCFieldSplit will attempt to automatically determine fields based on zero diagonal entries.
28947b752e3dSPatrick Sanan 
28957b752e3dSPatrick Sanan    Logically Collective
28967b752e3dSPatrick Sanan 
28977b752e3dSPatrick Sanan    Input Parameter:
28987b752e3dSPatrick Sanan .  pc   - the preconditioner context
28997b752e3dSPatrick Sanan 
29007b752e3dSPatrick Sanan    Output Parameter:
29017b752e3dSPatrick Sanan .  flg  - boolean indicating whether to detect fields or not
29027b752e3dSPatrick Sanan 
29037b752e3dSPatrick Sanan    Level: Intermediate
29047b752e3dSPatrick Sanan 
29057b752e3dSPatrick Sanan .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint()
29067b752e3dSPatrick Sanan 
29077b752e3dSPatrick Sanan @*/
29087b752e3dSPatrick Sanan PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC pc,PetscBool *flg)
29097b752e3dSPatrick Sanan {
29107b752e3dSPatrick Sanan   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
29117b752e3dSPatrick Sanan 
29127b752e3dSPatrick Sanan   PetscFunctionBegin;
29137b752e3dSPatrick Sanan   *flg = jac->detect;
29147b752e3dSPatrick Sanan   PetscFunctionReturn(0);
29157b752e3dSPatrick Sanan }
29167b752e3dSPatrick Sanan 
29177b752e3dSPatrick Sanan /*@
29187b752e3dSPatrick Sanan    PCFieldSplitSetDetectSaddlePoint - Sets flag indicating whether PCFieldSplit will attempt to automatically determine fields based on zero diagonal entries.
29197b752e3dSPatrick Sanan 
29207b752e3dSPatrick Sanan    Logically Collective
29217b752e3dSPatrick Sanan 
29227b752e3dSPatrick Sanan    Notes:
29237b752e3dSPatrick Sanan    Also sets the split type to PC_COMPOSITE_SCHUR (see PCFieldSplitSetType()) and the Schur preconditioner type to PC_FIELDSPLIT_SCHUR_PRE_SELF (see PCFieldSplitSetSchurPre()).
29247b752e3dSPatrick Sanan 
29257b752e3dSPatrick Sanan    Input Parameter:
29267b752e3dSPatrick Sanan .  pc   - the preconditioner context
29277b752e3dSPatrick Sanan 
29287b752e3dSPatrick Sanan    Output Parameter:
29297b752e3dSPatrick Sanan .  flg  - boolean indicating whether to detect fields or not
29307b752e3dSPatrick Sanan 
29317b752e3dSPatrick Sanan    Options Database Key:
29327b752e3dSPatrick Sanan .  -pc_fieldsplit_detect_saddle_point
29337b752e3dSPatrick Sanan 
29347b752e3dSPatrick Sanan    Level: Intermediate
29357b752e3dSPatrick Sanan 
29367b752e3dSPatrick Sanan .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint(), PCFieldSplitSetType(), PCFieldSplitSetSchurPre()
29377b752e3dSPatrick Sanan 
29387b752e3dSPatrick Sanan @*/
29397b752e3dSPatrick Sanan PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc,PetscBool flg)
29407b752e3dSPatrick Sanan {
29417b752e3dSPatrick Sanan   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
29427b752e3dSPatrick Sanan   PetscErrorCode ierr;
29437b752e3dSPatrick Sanan 
29447b752e3dSPatrick Sanan   PetscFunctionBegin;
29457b752e3dSPatrick Sanan   jac->detect = flg;
29467b752e3dSPatrick Sanan   if (jac->detect) {
29477b752e3dSPatrick Sanan     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
29487b752e3dSPatrick Sanan     ierr = PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELF,NULL);CHKERRQ(ierr);
29497b752e3dSPatrick Sanan   }
29507b752e3dSPatrick Sanan   PetscFunctionReturn(0);
29517b752e3dSPatrick Sanan }
29527b752e3dSPatrick Sanan 
29530971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
29540971522cSBarry Smith /*MC
2955a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2956a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
29570971522cSBarry Smith 
2958edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
2959edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
29600971522cSBarry Smith 
2961a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
296269a612a9SBarry Smith          and set the options directly on the resulting KSP object
29630971522cSBarry Smith 
29640971522cSBarry Smith    Level: intermediate
29650971522cSBarry Smith 
296679416396SBarry Smith    Options Database Keys:
296781540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
296881540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
296981540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
297081540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
2971a51937d4SCarola Kruse .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur,gkb> - type of relaxation or factorization splitting
2972a6a584a2SBarry Smith .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see PCFieldSplitSetSchurPre()
29737b752e3dSPatrick Sanan .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero diagonal and uses Schur complement with no preconditioner as the solver
297479416396SBarry Smith 
2975a51937d4SCarola Kruse .    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
29765d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
2977a51937d4SCarola Kruse -    Options prefix for inner solver when using Golub Kahan biadiagonalization preconditioner is -fieldsplit_0_
29785d4c12cdSJungho Lee 
2979c8a0d604SMatthew G Knepley    Notes:
2980c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
2981d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
2982d32f9abdSBarry Smith 
2983d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
2984d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
2985d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
2986d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
2987d32f9abdSBarry Smith 
2988c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
2989c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
2990c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
299113b05affSJed Brown $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2992c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
299313b05affSJed Brown      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
299413b05affSJed Brown $              S = A11 - A10 ksp(A00) A01
299513b05affSJed 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
2996686bed4dSStefano Zampini      in providing the SECOND split or 1 if not give). For PCFieldSplitGetSubKSP() when field number is 0,
2997c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2998a6a584a2SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.
29991d71051fSDmitry Karpeev 
3000a7476a74SDmitry Karpeev      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
30015668aaf4SBarry Smith      diag gives
3002c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
3003c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
3004c096484dSStefano 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
3005c096484dSStefano Zampini      can be turned off with PCFieldSplitSetSchurScale() or by command line -pc_fieldsplit_schur_scale 1.0. The lower factorization is the inverse of
3006c8a0d604SMatthew G Knepley $              (  A00   0 )
3007c8a0d604SMatthew G Knepley $              (  A10   S )
3008c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
3009c8a0d604SMatthew G Knepley $              ( A00 A01 )
3010c8a0d604SMatthew G Knepley $              (  0   S  )
3011c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
3012e69d4d44SBarry Smith 
3013edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
3014edf189efSBarry Smith      is used automatically for a second block.
3015edf189efSBarry Smith 
3016ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
3017ff218e97SBarry Smith      Generally it should be used with the AIJ format.
3018ff218e97SBarry Smith 
3019ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
3020ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
3021ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
30220716a85fSBarry Smith 
3023a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
30240971522cSBarry Smith 
30253bc631e0SBarry Smith    There is a nice discussion of block preconditioners in
30263bc631e0SBarry Smith 
30273bc631e0SBarry Smith [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
30283bc631e0SBarry Smith        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
30293bc631e0SBarry Smith        http://chess.cs.umd.edu/~elman/papers/tax.pdf
30303bc631e0SBarry Smith 
30314d308398SBarry Smith    The Constrained Pressure Preconditioner (CPR) can be implemented using PCCOMPOSITE with PCGALERKIN. CPR first solves an R A P subsystem, updates the
30324d308398SBarry Smith    residual on all variables (PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE)), and then applies a simple ILU like preconditioner on all the variables.
3033a6a584a2SBarry Smith 
3034a51937d4SCarola Kruse    The Golub-Kahan bidiagonalization preconditioner (gkb) can be applied to symmetric 2x2 block matrices of the shape
3035a51937d4SCarola Kruse $        ( A00  A01 )
3036a51937d4SCarola Kruse $        ( A01' 0   )
3037a51937d4SCarola Kruse    with A00 positive semi-definite. The implementation follows [Ar13]. Therein, we choose N := 1/nu * I and the (1,1)-block of the matrix is modified to H = A00 + nu*A01*A01'.
3038a51937d4SCarola Kruse    A linear system Hx = b has to be solved in each iteration of the GKB algorithm. The solver is chosen with the option prefix -fieldsplit_0_.
3039a51937d4SCarola Kruse 
3040a51937d4SCarola Kruse [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
3041a51937d4SCarola Kruse 
30427e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
3043285fb4e2SStefano Zampini            PCFieldSplitGetSubKSP(), PCFieldSplitSchurGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
30447b752e3dSPatrick Sanan           MatSchurComplementSetAinvType(), PCFieldSplitSetSchurScale(),
30457b752e3dSPatrick Sanan           PCFieldSplitSetDetectSaddlePoint()
30460971522cSBarry Smith M*/
30470971522cSBarry Smith 
30488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
30490971522cSBarry Smith {
30500971522cSBarry Smith   PetscErrorCode ierr;
30510971522cSBarry Smith   PC_FieldSplit  *jac;
30520971522cSBarry Smith 
30530971522cSBarry Smith   PetscFunctionBegin;
3054b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
30552fa5cd67SKarl Rupp 
30560971522cSBarry Smith   jac->bs                 = -1;
30570971522cSBarry Smith   jac->nsplits            = 0;
30583e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
3059e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
3060c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
3061c096484dSStefano Zampini   jac->schurscale         = -1.0;
3062fbe7908bSJed Brown   jac->dm_splits          = PETSC_TRUE;
30637b752e3dSPatrick Sanan   jac->detect             = PETSC_FALSE;
3064a51937d4SCarola Kruse   jac->gkbtol             = 1e-5;
3065a51937d4SCarola Kruse   jac->gkbdelay           = 5;
3066a51937d4SCarola Kruse   jac->gkbnu              = 1;
3067a51937d4SCarola Kruse   jac->gkbmaxit           = 100;
3068*de482cd7SCarola Kruse   jac->gkbmonitor         = PETSC_FALSE;
306951f519a2SBarry Smith 
30700971522cSBarry Smith   pc->data = (void*)jac;
30710971522cSBarry Smith 
30720971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
3073421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
30740971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
3075574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
30760971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
30770971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
30780971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
30790971522cSBarry Smith   pc->ops->applyrichardson = 0;
30800971522cSBarry Smith 
3081285fb4e2SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurGetSubKSP_C",PCFieldSplitSchurGetSubKSP_FieldSplit);CHKERRQ(ierr);
3082bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
3083bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
3084bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
3085bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
3086bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
30876dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);CHKERRQ(ierr);
30880971522cSBarry Smith   PetscFunctionReturn(0);
30890971522cSBarry Smith }
3090