xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h>  /*I "petscpc.h" I*/
2a80b646eSBarry Smith #include <petsc/private/kspimpl.h> /*  This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
31e25c274SJed Brown #include <petscdm.h>
40971522cSBarry Smith 
50a545947SLisandro Dalcin const char *const PCFieldSplitSchurPreTypes[]  = {"SELF", "SELFP", "A11", "USER", "FULL", "PCFieldSplitSchurPreType", "PC_FIELDSPLIT_SCHUR_PRE_", NULL};
60a545947SLisandro Dalcin const char *const PCFieldSplitSchurFactTypes[] = {"DIAG", "LOWER", "UPPER", "FULL", "PCFieldSplitSchurFactType", "PC_FIELDSPLIT_SCHUR_FACT_", NULL};
7c5d2311dSJed Brown 
89df09d43SBarry 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;
99df09d43SBarry Smith 
100971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
110971522cSBarry Smith struct _PC_FieldSplitLink {
1269a612a9SBarry Smith   KSP               ksp;
13443836d0SMatthew G Knepley   Vec               x, y, z;
14db4c96c1SJed Brown   char             *splitname;
150971522cSBarry Smith   PetscInt          nfields;
165d4c12cdSJungho Lee   PetscInt         *fields, *fields_col;
171b9fc7fcSBarry Smith   VecScatter        sctx;
18f5f0d762SBarry Smith   IS                is, is_col;
1951f519a2SBarry Smith   PC_FieldSplitLink next, previous;
209df09d43SBarry Smith   PetscLogEvent     event;
215ddf11f8SNicolas Barnafi 
225ddf11f8SNicolas Barnafi   /* Used only when setting coordinates with PCSetCoordinates */
235ddf11f8SNicolas Barnafi   PetscInt   dim;
245ddf11f8SNicolas Barnafi   PetscInt   ndofs;
255ddf11f8SNicolas Barnafi   PetscReal *coords;
260971522cSBarry Smith };
270971522cSBarry Smith 
280971522cSBarry Smith typedef struct {
2968dd23aaSBarry Smith   PCCompositeType type;
30ace3abfcSBarry Smith   PetscBool       defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
31ace3abfcSBarry Smith   PetscBool       splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
3230ad9308SMatthew Knepley   PetscInt        bs;           /* Block size for IS and Mat structures */
3330ad9308SMatthew Knepley   PetscInt        nsplits;      /* Number of field divisions defined */
3479416396SBarry Smith   Vec            *x, *y, w1, w2;
35519d70e2SJed Brown   Mat            *mat;    /* The diagonal block for each split */
36519d70e2SJed Brown   Mat            *pmat;   /* The preconditioning diagonal block for each split */
3730ad9308SMatthew Knepley   Mat            *Afield; /* The rows of the matrix associated with each split */
38ace3abfcSBarry Smith   PetscBool       issetup;
392fa5cd67SKarl Rupp 
4030ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4130ad9308SMatthew Knepley   Mat                       B;          /* The (0,1) block */
4230ad9308SMatthew Knepley   Mat                       C;          /* The (1,0) block */
43443836d0SMatthew G Knepley   Mat                       schur;      /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
44a7476a74SDmitry Karpeev   Mat                       schurp;     /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
45084e4875SJed Brown   Mat                       schur_user; /* User-provided preconditioning matrix for the Schur complement */
46084e4875SJed Brown   PCFieldSplitSchurPreType  schurpre;   /* Determines which preconditioning matrix is used for the Schur complement */
47c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
4830ad9308SMatthew Knepley   KSP                       kspschur;   /* The solver for S */
49443836d0SMatthew G Knepley   KSP                       kspupper;   /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
50c096484dSStefano Zampini   PetscScalar               schurscale; /* Scaling factor for the Schur complement solution with DIAG factorization */
51c096484dSStefano Zampini 
52a51937d4SCarola Kruse   /* Only used when Golub-Kahan bidiagonalization preconditioning is used */
53a51937d4SCarola Kruse   Mat          H;           /* The modified matrix H = A00 + nu*A01*A01'              */
54a51937d4SCarola Kruse   PetscReal    gkbtol;      /* Stopping tolerance for lower bound estimate            */
55a51937d4SCarola Kruse   PetscInt     gkbdelay;    /* The delay window for the stopping criterion            */
56a51937d4SCarola Kruse   PetscReal    gkbnu;       /* Parameter for augmented Lagrangian H = A + nu*A01*A01' */
57a51937d4SCarola Kruse   PetscInt     gkbmaxit;    /* Maximum number of iterations for outer loop            */
58de482cd7SCarola Kruse   PetscBool    gkbmonitor;  /* Monitor for gkb iterations and the lower bound error   */
59de482cd7SCarola Kruse   PetscViewer  gkbviewer;   /* Viewer context for gkbmonitor                          */
60e071a0a4SCarola Kruse   Vec          u, v, d, Hu; /* Work vectors for the GKB algorithm                     */
61e071a0a4SCarola Kruse   PetscScalar *vecz;        /* Contains intermediate values, eg for lower bound       */
62a51937d4SCarola Kruse 
6397bbdb24SBarry Smith   PC_FieldSplitLink head;
646dbb499eSCian Wilson   PetscBool         isrestrict;       /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
65c1570756SJed Brown   PetscBool         suboptionsset;    /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
664ab8060aSDmitry Karpeev   PetscBool         dm_splits;        /* Whether to use DMCreateFieldDecomposition() whenever possible */
67c84da90fSDmitry Karpeev   PetscBool         diag_use_amat;    /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
68c84da90fSDmitry Karpeev   PetscBool         offdiag_use_amat; /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
697b752e3dSPatrick Sanan   PetscBool         detect;           /* Whether to form 2-way split by finding zero diagonal entries */
705ddf11f8SNicolas Barnafi   PetscBool         coordinates_set;  /* Whether PCSetCoordinates has been called */
710971522cSBarry Smith } PC_FieldSplit;
720971522cSBarry Smith 
7316913363SBarry Smith /*
74f1580f4eSBarry Smith     Note:
7595452b02SPatrick Sanan     there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
7616913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
7716913363SBarry Smith    PC you could change this.
7816913363SBarry Smith */
79084e4875SJed Brown 
80e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
81084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
829371c9d4SSatish Balay static Mat FieldSplitSchurPre(PC_FieldSplit *jac) {
83084e4875SJed Brown   switch (jac->schurpre) {
84084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
85a7476a74SDmitry Karpeev   case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
86e87fab1bSBarry Smith   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
87e74569cdSMatthew G. Knepley   case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
88084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
899371c9d4SSatish Balay   default: return jac->schur_user ? jac->schur_user : jac->pmat[1];
90084e4875SJed Brown   }
91084e4875SJed Brown }
92084e4875SJed Brown 
939804daf3SBarry Smith #include <petscdraw.h>
949371c9d4SSatish Balay static PetscErrorCode PCView_FieldSplit(PC pc, PetscViewer viewer) {
950971522cSBarry Smith   PC_FieldSplit    *jac = (PC_FieldSplit *)pc->data;
96d9884438SBarry Smith   PetscBool         iascii, isdraw;
970971522cSBarry Smith   PetscInt          i, j;
985a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
990971522cSBarry Smith 
1000971522cSBarry Smith   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1030971522cSBarry Smith   if (iascii) {
1042c9966d7SBarry Smith     if (jac->bs > 0) {
10563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FieldSplit with %s composition: total splits = %" PetscInt_FMT ", blocksize = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits, jac->bs));
1062c9966d7SBarry Smith     } else {
10763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FieldSplit with %s composition: total splits = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits));
1082c9966d7SBarry Smith     }
10948a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for blocks\n"));
11048a46eb9SPierre Jolivet     if (jac->diag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for diagonal blocks\n"));
11148a46eb9SPierre Jolivet     if (jac->offdiag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for off-diagonal blocks\n"));
1129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Solver info for each split is in the following KSP objects:\n"));
1130971522cSBarry Smith     for (i = 0; i < jac->nsplits; i++) {
1141ab39975SBarry Smith       if (ilink->fields) {
11563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Fields ", i));
1169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1175a9f2f41SSatish Balay         for (j = 0; j < ilink->nfields; j++) {
11848a46eb9SPierre Jolivet           if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ","));
11963a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, ilink->fields[j]));
1200971522cSBarry Smith         }
1219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1231ab39975SBarry Smith       } else {
12463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Defined by IS\n", i));
1251ab39975SBarry Smith       }
1269566063dSJacob Faibussowitsch       PetscCall(KSPView(ilink->ksp, viewer));
1275a9f2f41SSatish Balay       ilink = ilink->next;
1280971522cSBarry Smith     }
1292fa5cd67SKarl Rupp   }
1302fa5cd67SKarl Rupp 
1312fa5cd67SKarl Rupp   if (isdraw) {
132d9884438SBarry Smith     PetscDraw draw;
133d9884438SBarry Smith     PetscReal x, y, w, wd;
134d9884438SBarry Smith 
1359566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1369566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
137d9884438SBarry Smith     w  = 2 * PetscMin(1.0 - x, x);
138d9884438SBarry Smith     wd = w / (jac->nsplits + 1);
139d9884438SBarry Smith     x  = x - wd * (jac->nsplits - 1) / 2.0;
140d9884438SBarry Smith     for (i = 0; i < jac->nsplits; i++) {
1419566063dSJacob Faibussowitsch       PetscCall(PetscDrawPushCurrentPoint(draw, x, y));
1429566063dSJacob Faibussowitsch       PetscCall(KSPView(ilink->ksp, viewer));
1439566063dSJacob Faibussowitsch       PetscCall(PetscDrawPopCurrentPoint(draw));
144d9884438SBarry Smith       x += wd;
145d9884438SBarry Smith       ilink = ilink->next;
146d9884438SBarry Smith     }
1470971522cSBarry Smith   }
1480971522cSBarry Smith   PetscFunctionReturn(0);
1490971522cSBarry Smith }
1500971522cSBarry Smith 
1519371c9d4SSatish Balay static PetscErrorCode PCView_FieldSplit_Schur(PC pc, PetscViewer viewer) {
1523b224e63SBarry Smith   PC_FieldSplit             *jac = (PC_FieldSplit *)pc->data;
1534996c5bdSBarry Smith   PetscBool                  iascii, isdraw;
1543b224e63SBarry Smith   PetscInt                   i, j;
1553b224e63SBarry Smith   PC_FieldSplitLink          ilink = jac->head;
156a9908d51SBarry Smith   MatSchurComplementAinvType atype;
1573b224e63SBarry Smith 
1583b224e63SBarry Smith   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1609566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1613b224e63SBarry Smith   if (iascii) {
1622c9966d7SBarry Smith     if (jac->bs > 0) {
16363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FieldSplit with Schur preconditioner, blocksize = %" PetscInt_FMT ", factorization %s\n", jac->bs, PCFieldSplitSchurFactTypes[jac->schurfactorization]));
1642c9966d7SBarry Smith     } else {
1659566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FieldSplit with Schur preconditioner, factorization %s\n", PCFieldSplitSchurFactTypes[jac->schurfactorization]));
1662c9966d7SBarry Smith     }
16748a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for blocks\n"));
1683e8b8b31SMatthew G Knepley     switch (jac->schurpre) {
1699371c9d4SSatish Balay     case PC_FIELDSPLIT_SCHUR_PRE_SELF: PetscCall(PetscViewerASCIIPrintf(viewer, "  Preconditioner for the Schur complement formed from S itself\n")); break;
170a7476a74SDmitry Karpeev     case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
1719566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementGetAinvType(jac->schur, &atype));
1729371c9d4SSatish Balay       PetscCall(PetscViewerASCIIPrintf(viewer, "  Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses A00's %sinverse\n", atype == MAT_SCHUR_COMPLEMENT_AINV_DIAG ? "diagonal's " : (atype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG ? "block diagonal's " : (atype == MAT_SCHUR_COMPLEMENT_AINV_FULL ? "full " : "lumped diagonal's "))));
173a9908d51SBarry Smith       break;
1749371c9d4SSatish Balay     case PC_FIELDSPLIT_SCHUR_PRE_A11: PetscCall(PetscViewerASCIIPrintf(viewer, "  Preconditioner for the Schur complement formed from A11\n")); break;
1759371c9d4SSatish Balay     case PC_FIELDSPLIT_SCHUR_PRE_FULL: PetscCall(PetscViewerASCIIPrintf(viewer, "  Preconditioner for the Schur complement formed from the exact Schur complement\n")); break;
1763e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1773e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1789566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Preconditioner for the Schur complement formed from user provided matrix\n"));
1793e8b8b31SMatthew G Knepley       } else {
1809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Preconditioner for the Schur complement formed from A11\n"));
1813e8b8b31SMatthew G Knepley       }
1823e8b8b31SMatthew G Knepley       break;
1839371c9d4SSatish Balay     default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1843e8b8b31SMatthew G Knepley     }
1859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Split info:\n"));
1869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1873b224e63SBarry Smith     for (i = 0; i < jac->nsplits; i++) {
1883b224e63SBarry Smith       if (ilink->fields) {
18963a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Fields ", i));
1909566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1913b224e63SBarry Smith         for (j = 0; j < ilink->nfields; j++) {
19248a46eb9SPierre Jolivet           if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ","));
19363a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, ilink->fields[j]));
1943b224e63SBarry Smith         }
1959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1969566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1973b224e63SBarry Smith       } else {
19863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Defined by IS\n", i));
1993b224e63SBarry Smith       }
2003b224e63SBarry Smith       ilink = ilink->next;
2013b224e63SBarry Smith     }
2029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for A00 block\n"));
2039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
20406de4afeSJed Brown     if (jac->head) {
2059566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->head->ksp, viewer));
2069566063dSJacob Faibussowitsch     } else PetscCall(PetscViewerASCIIPrintf(viewer, "  not yet available\n"));
2079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
20806de4afeSJed Brown     if (jac->head && jac->kspupper != jac->head->ksp) {
2099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for upper A00 in upper triangular factor \n"));
2109566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2119566063dSJacob Faibussowitsch       if (jac->kspupper) PetscCall(KSPView(jac->kspupper, viewer));
2129566063dSJacob Faibussowitsch       else PetscCall(PetscViewerASCIIPrintf(viewer, "  not yet available\n"));
2139566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
214443836d0SMatthew G Knepley     }
2159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for S = A11 - A10 inv(A00) A01 \n"));
2169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
21712cae6f2SJed Brown     if (jac->kspschur) {
2189566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->kspschur, viewer));
21912cae6f2SJed Brown     } else {
2209566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  not yet available\n"));
22112cae6f2SJed Brown     }
2229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
2239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
22406de4afeSJed Brown   } else if (isdraw && jac->head) {
2254996c5bdSBarry Smith     PetscDraw draw;
2264996c5bdSBarry Smith     PetscReal x, y, w, wd, h;
2274996c5bdSBarry Smith     PetscInt  cnt = 2;
2284996c5bdSBarry Smith     char      str[32];
2294996c5bdSBarry Smith 
2309566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2319566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
232c74581afSBarry Smith     if (jac->kspupper != jac->head->ksp) cnt++;
233c74581afSBarry Smith     w  = 2 * PetscMin(1.0 - x, x);
234c74581afSBarry Smith     wd = w / (cnt + 1);
235c74581afSBarry Smith 
2369566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, 32, "Schur fact. %s", PCFieldSplitSchurFactTypes[jac->schurfactorization]));
2379566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
2384996c5bdSBarry Smith     y -= h;
2394996c5bdSBarry Smith     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) {
2409566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 32, "Prec. for Schur from %s", PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]));
2413b224e63SBarry Smith     } else {
2429566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 32, "Prec. for Schur from %s", PCFieldSplitSchurPreTypes[jac->schurpre]));
2434996c5bdSBarry Smith     }
2449566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x + wd * (cnt - 1) / 2.0, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
2454996c5bdSBarry Smith     y -= h;
2464996c5bdSBarry Smith     x = x - wd * (cnt - 1) / 2.0;
2474996c5bdSBarry Smith 
2489566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, y));
2499566063dSJacob Faibussowitsch     PetscCall(KSPView(jac->head->ksp, viewer));
2509566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2514996c5bdSBarry Smith     if (jac->kspupper != jac->head->ksp) {
2524996c5bdSBarry Smith       x += wd;
2539566063dSJacob Faibussowitsch       PetscCall(PetscDrawPushCurrentPoint(draw, x, y));
2549566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->kspupper, viewer));
2559566063dSJacob Faibussowitsch       PetscCall(PetscDrawPopCurrentPoint(draw));
2564996c5bdSBarry Smith     }
2574996c5bdSBarry Smith     x += wd;
2589566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, y));
2599566063dSJacob Faibussowitsch     PetscCall(KSPView(jac->kspschur, viewer));
2609566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2613b224e63SBarry Smith   }
2623b224e63SBarry Smith   PetscFunctionReturn(0);
2633b224e63SBarry Smith }
2643b224e63SBarry Smith 
2659371c9d4SSatish Balay static PetscErrorCode PCView_FieldSplit_GKB(PC pc, PetscViewer viewer) {
266a51937d4SCarola Kruse   PC_FieldSplit    *jac = (PC_FieldSplit *)pc->data;
267a51937d4SCarola Kruse   PetscBool         iascii, isdraw;
268a51937d4SCarola Kruse   PetscInt          i, j;
269a51937d4SCarola Kruse   PC_FieldSplitLink ilink = jac->head;
270a51937d4SCarola Kruse 
271a51937d4SCarola Kruse   PetscFunctionBegin;
2729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
274a51937d4SCarola Kruse   if (iascii) {
275a51937d4SCarola Kruse     if (jac->bs > 0) {
27663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FieldSplit with %s composition: total splits = %" PetscInt_FMT ", blocksize = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits, jac->bs));
277a51937d4SCarola Kruse     } else {
27863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FieldSplit with %s composition: total splits = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits));
279a51937d4SCarola Kruse     }
28048a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for blocks\n"));
28148a46eb9SPierre Jolivet     if (jac->diag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for diagonal blocks\n"));
28248a46eb9SPierre Jolivet     if (jac->offdiag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat (not Pmat) as operator for off-diagonal blocks\n"));
283a51937d4SCarola Kruse 
28463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Stopping tolerance=%.1e, delay in error estimate=%" PetscInt_FMT ", maximum iterations=%" PetscInt_FMT "\n", (double)jac->gkbtol, jac->gkbdelay, jac->gkbmaxit));
2859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Solver info for H = A00 + nu*A01*A01' matrix:\n"));
2869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
287a51937d4SCarola Kruse 
288a51937d4SCarola Kruse     if (ilink->fields) {
28963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Split number 0 Fields "));
2909566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
291a51937d4SCarola Kruse       for (j = 0; j < ilink->nfields; j++) {
29248a46eb9SPierre Jolivet         if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ","));
29363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, ilink->fields[j]));
294a51937d4SCarola Kruse       }
2959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
2969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
297a51937d4SCarola Kruse     } else {
29863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Split number 0 Defined by IS\n"));
299a51937d4SCarola Kruse     }
3009566063dSJacob Faibussowitsch     PetscCall(KSPView(ilink->ksp, viewer));
301a51937d4SCarola Kruse 
3029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
303a51937d4SCarola Kruse   }
304a51937d4SCarola Kruse 
305a51937d4SCarola Kruse   if (isdraw) {
306a51937d4SCarola Kruse     PetscDraw draw;
307a51937d4SCarola Kruse     PetscReal x, y, w, wd;
308a51937d4SCarola Kruse 
3099566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
3109566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
311a51937d4SCarola Kruse     w  = 2 * PetscMin(1.0 - x, x);
312a51937d4SCarola Kruse     wd = w / (jac->nsplits + 1);
313a51937d4SCarola Kruse     x  = x - wd * (jac->nsplits - 1) / 2.0;
314a51937d4SCarola Kruse     for (i = 0; i < jac->nsplits; i++) {
3159566063dSJacob Faibussowitsch       PetscCall(PetscDrawPushCurrentPoint(draw, x, y));
3169566063dSJacob Faibussowitsch       PetscCall(KSPView(ilink->ksp, viewer));
3179566063dSJacob Faibussowitsch       PetscCall(PetscDrawPopCurrentPoint(draw));
318a51937d4SCarola Kruse       x += wd;
319a51937d4SCarola Kruse       ilink = ilink->next;
320a51937d4SCarola Kruse     }
321a51937d4SCarola Kruse   }
322a51937d4SCarola Kruse   PetscFunctionReturn(0);
323a51937d4SCarola Kruse }
324a51937d4SCarola Kruse 
3256c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
3269371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) {
3276c924f48SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
3285d4c12cdSJungho Lee   PetscInt       i, nfields, *ifields, nfields_col, *ifields_col;
3295d4c12cdSJungho Lee   PetscBool      flg, flg_col;
3305d4c12cdSJungho Lee   char           optionname[128], splitname[8], optionname_col[128];
3316c924f48SJed Brown 
3326c924f48SJed Brown   PetscFunctionBegin;
3339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jac->bs, &ifields));
3349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jac->bs, &ifields_col));
3356c924f48SJed Brown   for (i = 0, flg = PETSC_TRUE;; i++) {
33663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(splitname, sizeof(splitname), "%" PetscInt_FMT, i));
33763a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%" PetscInt_FMT "_fields", i));
33863a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(optionname_col, sizeof(optionname_col), "-pc_fieldsplit_%" PetscInt_FMT "_fields_col", i));
3396c924f48SJed Brown     nfields     = jac->bs;
34029499fbbSJungho Lee     nfields_col = jac->bs;
3419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetIntArray(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg));
3429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetIntArray(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, optionname_col, ifields_col, &nfields_col, &flg_col));
3436c924f48SJed Brown     if (!flg) break;
3445d4c12cdSJungho Lee     else if (flg && !flg_col) {
34528b400f6SJacob Faibussowitsch       PetscCheck(nfields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot list zero fields");
3469566063dSJacob Faibussowitsch       PetscCall(PCFieldSplitSetFields(pc, splitname, nfields, ifields, ifields));
3472fa5cd67SKarl Rupp     } else {
3487827d75bSBarry Smith       PetscCheck(nfields && nfields_col, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot list zero fields");
34908401ef6SPierre Jolivet       PetscCheck(nfields == nfields_col, PETSC_COMM_SELF, PETSC_ERR_USER, "Number of row and column fields must match");
3509566063dSJacob Faibussowitsch       PetscCall(PCFieldSplitSetFields(pc, splitname, nfields, ifields, ifields_col));
3515d4c12cdSJungho Lee     }
3526c924f48SJed Brown   }
3536c924f48SJed Brown   if (i > 0) {
3546c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
3556c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
3566c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
3576c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
3586c924f48SJed Brown   }
3599566063dSJacob Faibussowitsch   PetscCall(PetscFree(ifields));
3609566063dSJacob Faibussowitsch   PetscCall(PetscFree(ifields_col));
3616c924f48SJed Brown   PetscFunctionReturn(0);
3626c924f48SJed Brown }
3636c924f48SJed Brown 
3649371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetDefaults(PC pc) {
3650971522cSBarry Smith   PC_FieldSplit    *jac                = (PC_FieldSplit *)pc->data;
3665a9f2f41SSatish Balay   PC_FieldSplitLink ilink              = jac->head;
3677b752e3dSPatrick Sanan   PetscBool         fieldsplit_default = PETSC_FALSE, coupling = PETSC_FALSE;
3686c924f48SJed Brown   PetscInt          i;
3690971522cSBarry Smith 
3700971522cSBarry Smith   PetscFunctionBegin;
3717287d2a3SDmitry Karpeev   /*
372f5f0d762SBarry Smith    Kinda messy, but at least this now uses DMCreateFieldDecomposition().
3737287d2a3SDmitry Karpeev    Should probably be rewritten.
3747287d2a3SDmitry Karpeev    */
375f5f0d762SBarry Smith   if (!ilink) {
3769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_fieldsplit_detect_coupling", &coupling, NULL));
3777b752e3dSPatrick Sanan     if (pc->dm && jac->dm_splits && !jac->detect && !coupling) {
378bafc1b83SMatthew G Knepley       PetscInt  numFields, f, i, j;
3790784a22cSJed Brown       char    **fieldNames;
3807b62db95SJungho Lee       IS       *fields;
381e7c4fc90SDmitry Karpeev       DM       *dms;
382bafc1b83SMatthew G Knepley       DM        subdm[128];
383bafc1b83SMatthew G Knepley       PetscBool flg;
384bafc1b83SMatthew G Knepley 
3859566063dSJacob Faibussowitsch       PetscCall(DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms));
386bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
387bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE;; i++) {
388bafc1b83SMatthew G Knepley         PetscInt ifields[128];
389bafc1b83SMatthew G Knepley         IS       compField;
390bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
391bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
392bafc1b83SMatthew G Knepley 
39363a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%" PetscInt_FMT "_fields", i));
3949566063dSJacob Faibussowitsch         PetscCall(PetscOptionsGetIntArray(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg));
395bafc1b83SMatthew G Knepley         if (!flg) break;
39663a3b9bcSJacob Faibussowitsch         PetscCheck(numFields <= 128, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot currently support %" PetscInt_FMT " > 128 fields", numFields);
3979566063dSJacob Faibussowitsch         PetscCall(DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]));
398bafc1b83SMatthew G Knepley         if (nfields == 1) {
3999566063dSJacob Faibussowitsch           PetscCall(PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField));
400bafc1b83SMatthew G Knepley         } else {
40163a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(splitname, sizeof(splitname), "%" PetscInt_FMT, i));
4029566063dSJacob Faibussowitsch           PetscCall(PCFieldSplitSetIS(pc, splitname, compField));
4037287d2a3SDmitry Karpeev         }
4049566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&compField));
405bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
406bafc1b83SMatthew G Knepley           f = ifields[j];
4079566063dSJacob Faibussowitsch           PetscCall(PetscFree(fieldNames[f]));
4089566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&fields[f]));
4097b62db95SJungho Lee         }
410bafc1b83SMatthew G Knepley       }
411bafc1b83SMatthew G Knepley       if (i == 0) {
412bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
4139566063dSJacob Faibussowitsch           PetscCall(PCFieldSplitSetIS(pc, fieldNames[f], fields[f]));
4149566063dSJacob Faibussowitsch           PetscCall(PetscFree(fieldNames[f]));
4159566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&fields[f]));
416bafc1b83SMatthew G Knepley         }
417bafc1b83SMatthew G Knepley       } else {
41848a46eb9SPierre Jolivet         for (j = 0; j < numFields; j++) PetscCall(DMDestroy(dms + j));
4199566063dSJacob Faibussowitsch         PetscCall(PetscFree(dms));
4209566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(i, &dms));
4212fa5cd67SKarl Rupp         for (j = 0; j < i; ++j) dms[j] = subdm[j];
422bafc1b83SMatthew G Knepley       }
4239566063dSJacob Faibussowitsch       PetscCall(PetscFree(fieldNames));
4249566063dSJacob Faibussowitsch       PetscCall(PetscFree(fields));
425e7c4fc90SDmitry Karpeev       if (dms) {
4269566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n"));
427bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
4287287d2a3SDmitry Karpeev           const char *prefix;
4299566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp), &prefix));
4309566063dSJacob Faibussowitsch           PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix));
4319566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(ilink->ksp, dms[i]));
4329566063dSJacob Faibussowitsch           PetscCall(KSPSetDMActive(ilink->ksp, PETSC_FALSE));
43395dbaa6fSMatthew G. Knepley           {
43495dbaa6fSMatthew G. Knepley             PetscErrorCode (*func)(KSP, Mat, Mat, void *);
43595dbaa6fSMatthew G. Knepley             void *ctx;
43695dbaa6fSMatthew G. Knepley 
4379566063dSJacob Faibussowitsch             PetscCall(DMKSPGetComputeOperators(pc->dm, &func, &ctx));
4389566063dSJacob Faibussowitsch             PetscCall(DMKSPSetComputeOperators(dms[i], func, ctx));
43995dbaa6fSMatthew G. Knepley           }
4409566063dSJacob Faibussowitsch           PetscCall(PetscObjectIncrementTabLevel((PetscObject)dms[i], (PetscObject)ilink->ksp, 0));
4419566063dSJacob Faibussowitsch           PetscCall(DMDestroy(&dms[i]));
4422fa5ba8aSJed Brown         }
4439566063dSJacob Faibussowitsch         PetscCall(PetscFree(dms));
4448b8307b2SJed Brown       }
44566ffff09SJed Brown     } else {
446521d7252SBarry Smith       if (jac->bs <= 0) {
447704ba839SBarry Smith         if (pc->pmat) {
4489566063dSJacob Faibussowitsch           PetscCall(MatGetBlockSize(pc->pmat, &jac->bs));
4492fa5cd67SKarl Rupp         } else jac->bs = 1;
450521d7252SBarry Smith       }
451d32f9abdSBarry Smith 
4527b752e3dSPatrick Sanan       if (jac->detect) {
4536ce1633cSBarry Smith         IS       zerodiags, rest;
4546ce1633cSBarry Smith         PetscInt nmin, nmax;
4556ce1633cSBarry Smith 
4569566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(pc->mat, &nmin, &nmax));
4577199da05SBarry Smith         if (jac->diag_use_amat) {
4589566063dSJacob Faibussowitsch           PetscCall(MatFindZeroDiagonals(pc->mat, &zerodiags));
4597199da05SBarry Smith         } else {
4609566063dSJacob Faibussowitsch           PetscCall(MatFindZeroDiagonals(pc->pmat, &zerodiags));
4617199da05SBarry Smith         }
4629566063dSJacob Faibussowitsch         PetscCall(ISComplement(zerodiags, nmin, nmax, &rest));
4639566063dSJacob Faibussowitsch         PetscCall(PCFieldSplitSetIS(pc, "0", rest));
4649566063dSJacob Faibussowitsch         PetscCall(PCFieldSplitSetIS(pc, "1", zerodiags));
4659566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&zerodiags));
4669566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&rest));
4673a062f41SBarry Smith       } else if (coupling) {
4683a062f41SBarry Smith         IS       coupling, rest;
4693a062f41SBarry Smith         PetscInt nmin, nmax;
4703a062f41SBarry Smith 
4719566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(pc->mat, &nmin, &nmax));
4727199da05SBarry Smith         if (jac->offdiag_use_amat) {
4739566063dSJacob Faibussowitsch           PetscCall(MatFindOffBlockDiagonalEntries(pc->mat, &coupling));
4747199da05SBarry Smith         } else {
4759566063dSJacob Faibussowitsch           PetscCall(MatFindOffBlockDiagonalEntries(pc->pmat, &coupling));
4767199da05SBarry Smith         }
4779566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc->mat), nmax - nmin, nmin, 1, &rest));
4789566063dSJacob Faibussowitsch         PetscCall(ISSetIdentity(rest));
4799566063dSJacob Faibussowitsch         PetscCall(PCFieldSplitSetIS(pc, "0", rest));
4809566063dSJacob Faibussowitsch         PetscCall(PCFieldSplitSetIS(pc, "1", coupling));
4819566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&coupling));
4829566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&rest));
4836ce1633cSBarry Smith       } else {
4849566063dSJacob Faibussowitsch         PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_fieldsplit_default", &fieldsplit_default, NULL));
4857287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
486d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
487d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
4889566063dSJacob Faibussowitsch           PetscCall(PCFieldSplitSetRuntimeSplits_Private(pc));
4899566063dSJacob Faibussowitsch           if (jac->splitdefined) PetscCall(PetscInfo(pc, "Splits defined using the options database\n"));
490d32f9abdSBarry Smith         }
4916dbb499eSCian Wilson         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
4929f001fe8SStefano Zampini           Mat       M = pc->pmat;
493f3b928b9SStefano Zampini           PetscBool isnest;
494f3b928b9SStefano Zampini 
4959566063dSJacob Faibussowitsch           PetscCall(PetscInfo(pc, "Using default splitting of fields\n"));
4969566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATNEST, &isnest));
497f3b928b9SStefano Zampini           if (!isnest) {
4989f001fe8SStefano Zampini             M = pc->mat;
4999566063dSJacob Faibussowitsch             PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat, MATNEST, &isnest));
500f3b928b9SStefano Zampini           }
501f3b928b9SStefano Zampini           if (isnest) {
502f3b928b9SStefano Zampini             IS      *fields;
503f3b928b9SStefano Zampini             PetscInt nf;
504f3b928b9SStefano Zampini 
5059566063dSJacob Faibussowitsch             PetscCall(MatNestGetSize(M, &nf, NULL));
5069566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(nf, &fields));
5079566063dSJacob Faibussowitsch             PetscCall(MatNestGetISs(M, fields, NULL));
50848a46eb9SPierre Jolivet             for (i = 0; i < nf; i++) PetscCall(PCFieldSplitSetIS(pc, NULL, fields[i]));
5099566063dSJacob Faibussowitsch             PetscCall(PetscFree(fields));
510f3b928b9SStefano Zampini           } else {
511db4c96c1SJed Brown             for (i = 0; i < jac->bs; i++) {
5126c924f48SJed Brown               char splitname[8];
51363a3b9bcSJacob Faibussowitsch               PetscCall(PetscSNPrintf(splitname, sizeof(splitname), "%" PetscInt_FMT, i));
5149566063dSJacob Faibussowitsch               PetscCall(PCFieldSplitSetFields(pc, splitname, 1, &i, &i));
51579416396SBarry Smith             }
5165d4c12cdSJungho Lee             jac->defaultsplit = PETSC_TRUE;
517521d7252SBarry Smith           }
51866ffff09SJed Brown         }
5196ce1633cSBarry Smith       }
520f3b928b9SStefano Zampini     }
521edf189efSBarry Smith   } else if (jac->nsplits == 1) {
522edf189efSBarry Smith     if (ilink->is) {
523edf189efSBarry Smith       IS       is2;
524edf189efSBarry Smith       PetscInt nmin, nmax;
525edf189efSBarry Smith 
5269566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(pc->mat, &nmin, &nmax));
5279566063dSJacob Faibussowitsch       PetscCall(ISComplement(ilink->is, nmin, nmax, &is2));
5289566063dSJacob Faibussowitsch       PetscCall(PCFieldSplitSetIS(pc, "1", is2));
5299566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is2));
53082f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Must provide at least two sets of fields to PCFieldSplit()");
531edf189efSBarry Smith   }
532d0af7cd3SBarry Smith 
53363a3b9bcSJacob Faibussowitsch   PetscCheck(jac->nsplits >= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Unhandled case, must have at least two fields, not %" PetscInt_FMT, jac->nsplits);
53469a612a9SBarry Smith   PetscFunctionReturn(0);
53569a612a9SBarry Smith }
53669a612a9SBarry Smith 
5379371c9d4SSatish Balay static PetscErrorCode MatGolubKahanComputeExplicitOperator(Mat A, Mat B, Mat C, Mat *H, PetscReal gkbnu) {
538a51937d4SCarola Kruse   Mat       BT, T;
539de482cd7SCarola Kruse   PetscReal nrmT, nrmB;
540a51937d4SCarola Kruse 
541a51937d4SCarola Kruse   PetscFunctionBegin;
5429566063dSJacob Faibussowitsch   PetscCall(MatHermitianTranspose(C, MAT_INITIAL_MATRIX, &T)); /* Test if augmented matrix is symmetric */
5439566063dSJacob Faibussowitsch   PetscCall(MatAXPY(T, -1.0, B, DIFFERENT_NONZERO_PATTERN));
5449566063dSJacob Faibussowitsch   PetscCall(MatNorm(T, NORM_1, &nrmT));
5459566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_1, &nrmB));
546049d1499SBarry Smith   PetscCheck(nrmB <= 0 || nrmT / nrmB < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix is not symmetric/hermitian, GKB is not applicable.");
547049d1499SBarry Smith 
548a51937d4SCarola Kruse   /* Compute augmented Lagrangian matrix H = A00 + nu*A01*A01'. This corresponds to */
549a51937d4SCarola Kruse   /* setting N := 1/nu*I in [Ar13].                                                 */
5509566063dSJacob Faibussowitsch   PetscCall(MatHermitianTranspose(B, MAT_INITIAL_MATRIX, &BT));
5519566063dSJacob Faibussowitsch   PetscCall(MatMatMult(B, BT, MAT_INITIAL_MATRIX, PETSC_DEFAULT, H)); /* H = A01*A01'          */
5529566063dSJacob Faibussowitsch   PetscCall(MatAYPX(*H, gkbnu, A, DIFFERENT_NONZERO_PATTERN));        /* H = A00 + nu*A01*A01' */
553a51937d4SCarola Kruse 
5549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BT));
5559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&T));
556a51937d4SCarola Kruse   PetscFunctionReturn(0);
557a51937d4SCarola Kruse }
558a51937d4SCarola Kruse 
5592d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions, const char pre[], const char name[], const char *value[], PetscBool *flg);
560514bf10dSMatthew G Knepley 
5619371c9d4SSatish Balay static PetscErrorCode PCSetUp_FieldSplit(PC pc) {
56269a612a9SBarry Smith   PC_FieldSplit    *jac = (PC_FieldSplit *)pc->data;
5635a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
5642c9966d7SBarry Smith   PetscInt          i, nsplit;
5655f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
56669a612a9SBarry Smith 
56769a612a9SBarry Smith   PetscFunctionBegin;
5685da88fe4STristan Konolige   pc->failedreason = PC_NOERROR;
5699566063dSJacob Faibussowitsch   PetscCall(PCFieldSplitSetDefaults(pc));
57097bbdb24SBarry Smith   nsplit = jac->nsplits;
5715a9f2f41SSatish Balay   ilink  = jac->head;
57297bbdb24SBarry Smith 
57397bbdb24SBarry Smith   /* get the matrices for each split */
574704ba839SBarry Smith   if (!jac->issetup) {
5751b9fc7fcSBarry Smith     PetscInt rstart, rend, nslots, bs;
57697bbdb24SBarry Smith 
577704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
578704ba839SBarry Smith 
5795d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
5802c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
5812c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
5822c9966d7SBarry Smith     }
58351f519a2SBarry Smith     bs = jac->bs;
5849566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(pc->pmat, &rstart, &rend));
5851b9fc7fcSBarry Smith     nslots = (rend - rstart) / bs;
5861b9fc7fcSBarry Smith     for (i = 0; i < nsplit; i++) {
5871b9fc7fcSBarry Smith       if (jac->defaultsplit) {
5889566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), nslots, rstart + i, nsplit, &ilink->is));
5899566063dSJacob Faibussowitsch         PetscCall(ISDuplicate(ilink->is, &ilink->is_col));
590704ba839SBarry Smith       } else if (!ilink->is) {
591ccb205f8SBarry Smith         if (ilink->nfields > 1) {
5925f4ab4e1SJungho Lee           PetscInt *ii, *jj, j, k, nfields = ilink->nfields, *fields = ilink->fields, *fields_col = ilink->fields_col;
5939566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(ilink->nfields * nslots, &ii));
5949566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(ilink->nfields * nslots, &jj));
5951b9fc7fcSBarry Smith           for (j = 0; j < nslots; j++) {
5961b9fc7fcSBarry Smith             for (k = 0; k < nfields; k++) {
5971b9fc7fcSBarry Smith               ii[nfields * j + k] = rstart + bs * j + fields[k];
5985f4ab4e1SJungho Lee               jj[nfields * j + k] = rstart + bs * j + fields_col[k];
59997bbdb24SBarry Smith             }
60097bbdb24SBarry Smith           }
6019566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), nslots * nfields, ii, PETSC_OWN_POINTER, &ilink->is));
6029566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), nslots * nfields, jj, PETSC_OWN_POINTER, &ilink->is_col));
6039566063dSJacob Faibussowitsch           PetscCall(ISSetBlockSize(ilink->is, nfields));
6049566063dSJacob Faibussowitsch           PetscCall(ISSetBlockSize(ilink->is_col, nfields));
605ccb205f8SBarry Smith         } else {
6069566063dSJacob Faibussowitsch           PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), nslots, rstart + ilink->fields[0], bs, &ilink->is));
6079566063dSJacob Faibussowitsch           PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), nslots, rstart + ilink->fields_col[0], bs, &ilink->is_col));
608ccb205f8SBarry Smith         }
6093e197d65SBarry Smith       }
6109566063dSJacob Faibussowitsch       PetscCall(ISSorted(ilink->is, &sorted));
6119566063dSJacob Faibussowitsch       if (ilink->is_col) PetscCall(ISSorted(ilink->is_col, &sorted_col));
6127827d75bSBarry Smith       PetscCheck(sorted && sorted_col, PETSC_COMM_SELF, PETSC_ERR_USER, "Fields must be sorted when creating split");
613704ba839SBarry Smith       ilink = ilink->next;
6141b9fc7fcSBarry Smith     }
6151b9fc7fcSBarry Smith   }
6161b9fc7fcSBarry Smith 
617704ba839SBarry Smith   ilink = jac->head;
61897bbdb24SBarry Smith   if (!jac->pmat) {
619bdddcaaaSMatthew G Knepley     Vec xtmp;
620bdddcaaaSMatthew G Knepley 
6219566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &xtmp, NULL));
6229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsplit, &jac->pmat));
6239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nsplit, &jac->x, nsplit, &jac->y));
624cf502942SBarry Smith     for (i = 0; i < nsplit; i++) {
625bdddcaaaSMatthew G Knepley       MatNullSpace sp;
626bdddcaaaSMatthew G Knepley 
627a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
6289566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject)ilink->is, "pmat", (PetscObject *)&jac->pmat[i]));
629a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
6309566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)jac->pmat[i]));
631a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
632a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
6332fa5cd67SKarl Rupp 
6349566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)jac->schur_user));
635a3df900dSMatthew G Knepley         }
636a3df900dSMatthew G Knepley       } else {
6373a062f41SBarry Smith         const char *prefix;
6389566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ilink->is_col, MAT_INITIAL_MATRIX, &jac->pmat[i]));
6399566063dSJacob Faibussowitsch         PetscCall(KSPGetOptionsPrefix(ilink->ksp, &prefix));
6409566063dSJacob Faibussowitsch         PetscCall(MatSetOptionsPrefix(jac->pmat[i], prefix));
6419566063dSJacob Faibussowitsch         PetscCall(MatViewFromOptions(jac->pmat[i], NULL, "-mat_view"));
642a3df900dSMatthew G Knepley       }
643bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
6449566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(jac->pmat[i], &jac->x[i], &jac->y[i]));
6459371c9d4SSatish Balay       ilink->x = jac->x[i];
6469371c9d4SSatish Balay       ilink->y = jac->y[i];
6479371c9d4SSatish Balay       ilink->z = NULL;
648bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
6499566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(xtmp, ilink->is, jac->x[i], NULL, &ilink->sctx));
6509566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject)ilink->is, "nearnullspace", (PetscObject *)&sp));
65148a46eb9SPierre Jolivet       if (sp) PetscCall(MatSetNearNullSpace(jac->pmat[i], sp));
652704ba839SBarry Smith       ilink = ilink->next;
653cf502942SBarry Smith     }
6549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xtmp));
65597bbdb24SBarry Smith   } else {
656ef7efd37SHong Zhang     MatReuse scall;
657ef7efd37SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
65848a46eb9SPierre Jolivet       for (i = 0; i < nsplit; i++) PetscCall(MatDestroy(&jac->pmat[i]));
659ef7efd37SHong Zhang       scall = MAT_INITIAL_MATRIX;
660ef7efd37SHong Zhang     } else scall = MAT_REUSE_MATRIX;
661ef7efd37SHong Zhang 
662cf502942SBarry Smith     for (i = 0; i < nsplit; i++) {
663a3df900dSMatthew G Knepley       Mat pmat;
664a3df900dSMatthew G Knepley 
665a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
6669566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject)ilink->is, "pmat", (PetscObject *)&pmat));
66748a46eb9SPierre Jolivet       if (!pmat) PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ilink->is_col, scall, &jac->pmat[i]));
668704ba839SBarry Smith       ilink = ilink->next;
669cf502942SBarry Smith     }
67097bbdb24SBarry Smith   }
6714e39094bSDmitry Karpeev   if (jac->diag_use_amat) {
672519d70e2SJed Brown     ilink = jac->head;
673519d70e2SJed Brown     if (!jac->mat) {
6749566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nsplit, &jac->mat));
675519d70e2SJed Brown       for (i = 0; i < nsplit; i++) {
6769566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ilink->is_col, MAT_INITIAL_MATRIX, &jac->mat[i]));
677519d70e2SJed Brown         ilink = ilink->next;
678519d70e2SJed Brown       }
679519d70e2SJed Brown     } else {
680ef7efd37SHong Zhang       MatReuse scall;
681ef7efd37SHong Zhang       if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
68248a46eb9SPierre Jolivet         for (i = 0; i < nsplit; i++) PetscCall(MatDestroy(&jac->mat[i]));
683ef7efd37SHong Zhang         scall = MAT_INITIAL_MATRIX;
684ef7efd37SHong Zhang       } else scall = MAT_REUSE_MATRIX;
685ef7efd37SHong Zhang 
686ef7efd37SHong Zhang       for (i = 0; i < nsplit; i++) {
6879566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ilink->is_col, scall, &jac->mat[i]));
688519d70e2SJed Brown         ilink = ilink->next;
689519d70e2SJed Brown       }
690519d70e2SJed Brown     }
691519d70e2SJed Brown   } else {
692519d70e2SJed Brown     jac->mat = jac->pmat;
693519d70e2SJed Brown   }
69497bbdb24SBarry Smith 
69553935eafSBarry Smith   /* Check for null space attached to IS */
69653935eafSBarry Smith   ilink = jac->head;
69753935eafSBarry Smith   for (i = 0; i < nsplit; i++) {
69853935eafSBarry Smith     MatNullSpace sp;
69953935eafSBarry Smith 
7009566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)ilink->is, "nullspace", (PetscObject *)&sp));
70148a46eb9SPierre Jolivet     if (sp) PetscCall(MatSetNullSpace(jac->mat[i], sp));
70253935eafSBarry Smith     ilink = ilink->next;
70353935eafSBarry Smith   }
70453935eafSBarry Smith 
705a51937d4SCarola Kruse   if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR && jac->type != PC_COMPOSITE_GKB) {
70668dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
7074e39094bSDmitry Karpeev     /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
70868dd23aaSBarry Smith     ilink = jac->head;
709e52d2c62SBarry Smith     if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
710e52d2c62SBarry 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 */
711e52d2c62SBarry Smith       if (!jac->Afield) {
7129566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(nsplit, &jac->Afield));
71380c96bb1SFande Kong         if (jac->offdiag_use_amat) {
7149566063dSJacob Faibussowitsch           PetscCall(MatCreateSubMatrix(pc->mat, ilink->next->is, ilink->is, MAT_INITIAL_MATRIX, &jac->Afield[1]));
715e52d2c62SBarry Smith         } else {
7169566063dSJacob Faibussowitsch           PetscCall(MatCreateSubMatrix(pc->pmat, ilink->next->is, ilink->is, MAT_INITIAL_MATRIX, &jac->Afield[1]));
71780c96bb1SFande Kong         }
71880c96bb1SFande Kong       } else {
719ef7efd37SHong Zhang         MatReuse scall;
720e9422dd5SStefano Zampini 
721ef7efd37SHong Zhang         if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
7229566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&jac->Afield[1]));
723ef7efd37SHong Zhang           scall = MAT_INITIAL_MATRIX;
724ef7efd37SHong Zhang         } else scall = MAT_REUSE_MATRIX;
725ef7efd37SHong Zhang 
72680c96bb1SFande Kong         if (jac->offdiag_use_amat) {
7279566063dSJacob Faibussowitsch           PetscCall(MatCreateSubMatrix(pc->mat, ilink->next->is, ilink->is, scall, &jac->Afield[1]));
72880c96bb1SFande Kong         } else {
7299566063dSJacob Faibussowitsch           PetscCall(MatCreateSubMatrix(pc->pmat, ilink->next->is, ilink->is, scall, &jac->Afield[1]));
73080c96bb1SFande Kong         }
731e52d2c62SBarry Smith       }
732e52d2c62SBarry Smith     } else {
73368dd23aaSBarry Smith       if (!jac->Afield) {
7349566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nsplit, &jac->Afield));
73568dd23aaSBarry Smith         for (i = 0; i < nsplit; i++) {
73680c96bb1SFande Kong           if (jac->offdiag_use_amat) {
7379566063dSJacob Faibussowitsch             PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, NULL, MAT_INITIAL_MATRIX, &jac->Afield[i]));
73880c96bb1SFande Kong           } else {
7399566063dSJacob Faibussowitsch             PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, NULL, MAT_INITIAL_MATRIX, &jac->Afield[i]));
74080c96bb1SFande Kong           }
74168dd23aaSBarry Smith           ilink = ilink->next;
74268dd23aaSBarry Smith         }
74368dd23aaSBarry Smith       } else {
744ef7efd37SHong Zhang         MatReuse scall;
745ef7efd37SHong Zhang         if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
74648a46eb9SPierre Jolivet           for (i = 0; i < nsplit; i++) PetscCall(MatDestroy(&jac->Afield[i]));
747ef7efd37SHong Zhang           scall = MAT_INITIAL_MATRIX;
748ef7efd37SHong Zhang         } else scall = MAT_REUSE_MATRIX;
749ef7efd37SHong Zhang 
75068dd23aaSBarry Smith         for (i = 0; i < nsplit; i++) {
75180c96bb1SFande Kong           if (jac->offdiag_use_amat) {
7529566063dSJacob Faibussowitsch             PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, NULL, scall, &jac->Afield[i]));
75380c96bb1SFande Kong           } else {
7549566063dSJacob Faibussowitsch             PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, NULL, scall, &jac->Afield[i]));
75580c96bb1SFande Kong           }
75668dd23aaSBarry Smith           ilink = ilink->next;
75768dd23aaSBarry Smith         }
75868dd23aaSBarry Smith       }
75968dd23aaSBarry Smith     }
760e52d2c62SBarry Smith   }
76168dd23aaSBarry Smith 
7623b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
7633b224e63SBarry Smith     IS          ccis;
764b94d7dedSBarry Smith     PetscBool   isset, isspd;
7654aa3045dSJed Brown     PetscInt    rstart, rend;
766093c86ffSJed Brown     char        lscname[256];
767093c86ffSJed Brown     PetscObject LSC_L;
768ce94432eSBarry Smith 
76908401ef6SPierre Jolivet     PetscCheck(nsplit == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "To use Schur complement preconditioner you must have exactly 2 fields");
77068dd23aaSBarry Smith 
771c096484dSStefano Zampini     /* If pc->mat is SPD, don't scale by -1 the Schur complement */
772c096484dSStefano Zampini     if (jac->schurscale == (PetscScalar)-1.0) {
773b94d7dedSBarry Smith       PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
774b94d7dedSBarry Smith       jac->schurscale = (isset && isspd) ? 1.0 : -1.0;
775c096484dSStefano Zampini     }
776c096484dSStefano Zampini 
777e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
7789566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(pc->mat, &rstart, &rend));
779e6cab6aaSJed Brown 
7803b224e63SBarry Smith     if (jac->schur) {
7810298fd71SBarry Smith       KSP      kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
782e9422dd5SStefano Zampini       MatReuse scall;
783e9422dd5SStefano Zampini 
784e9422dd5SStefano Zampini       if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
785e9422dd5SStefano Zampini         scall = MAT_INITIAL_MATRIX;
7869566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&jac->B));
7879566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&jac->C));
788e9422dd5SStefano Zampini       } else scall = MAT_REUSE_MATRIX;
789443836d0SMatthew G Knepley 
7909566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementGetKSP(jac->schur, &kspInner));
7913b224e63SBarry Smith       ilink = jac->head;
7929566063dSJacob Faibussowitsch       PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis));
7934e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
7949566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, scall, &jac->B));
795c84da90fSDmitry Karpeev       } else {
7969566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, scall, &jac->B));
797c84da90fSDmitry Karpeev       }
7989566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&ccis));
7993b224e63SBarry Smith       ilink = ilink->next;
8009566063dSJacob Faibussowitsch       PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis));
8014e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8029566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, scall, &jac->C));
803c84da90fSDmitry Karpeev       } else {
8049566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, scall, &jac->C));
805c84da90fSDmitry Karpeev       }
8069566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&ccis));
8079566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementUpdateSubMatrices(jac->schur, jac->mat[0], jac->pmat[0], jac->B, jac->C, jac->mat[1]));
808a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
8099566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&jac->schurp));
8109566063dSJacob Faibussowitsch         PetscCall(MatSchurComplementGetPmat(jac->schur, MAT_INITIAL_MATRIX, &jac->schurp));
811a7476a74SDmitry Karpeev       }
81248a46eb9SPierre Jolivet       if (kspA != kspInner) PetscCall(KSPSetOperators(kspA, jac->mat[0], jac->pmat[0]));
81348a46eb9SPierre Jolivet       if (kspUpper != kspA) PetscCall(KSPSetOperators(kspUpper, jac->mat[0], jac->pmat[0]));
8149566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->kspschur, jac->schur, FieldSplitSchurPre(jac)));
8153b224e63SBarry Smith     } else {
816bafc1b83SMatthew G Knepley       const char  *Dprefix;
817470b340bSDmitry Karpeev       char         schurprefix[256], schurmatprefix[256];
818514bf10dSMatthew G Knepley       char         schurtestoption[256];
819bdddcaaaSMatthew G Knepley       MatNullSpace sp;
820514bf10dSMatthew G Knepley       PetscBool    flg;
821686bed4dSStefano Zampini       KSP          kspt;
8223b224e63SBarry Smith 
823a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
8243b224e63SBarry Smith       ilink = jac->head;
8259566063dSJacob Faibussowitsch       PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis));
8264e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8279566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B));
828c84da90fSDmitry Karpeev       } else {
8299566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B));
830c84da90fSDmitry Karpeev       }
8319566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&ccis));
8323b224e63SBarry Smith       ilink = ilink->next;
8339566063dSJacob Faibussowitsch       PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis));
8344e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
8359566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C));
836c84da90fSDmitry Karpeev       } else {
8379566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C));
838c84da90fSDmitry Karpeev       }
8399566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&ccis));
84020252d06SBarry Smith 
841f5236f50SJed Brown       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
8429566063dSJacob Faibussowitsch       PetscCall(MatCreate(((PetscObject)jac->mat[0])->comm, &jac->schur));
8439566063dSJacob Faibussowitsch       PetscCall(MatSetType(jac->schur, MATSCHURCOMPLEMENT));
8449566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementSetSubMatrices(jac->schur, jac->mat[0], jac->pmat[0], jac->B, jac->C, jac->mat[1]));
8459566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname));
8469566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(jac->schur, schurmatprefix));
8479566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementGetKSP(jac->schur, &kspt));
8489566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(kspt, schurmatprefix));
849686bed4dSStefano Zampini 
850686bed4dSStefano Zampini       /* Note: this is not true in general */
8519566063dSJacob Faibussowitsch       PetscCall(MatGetNullSpace(jac->mat[1], &sp));
8521baa6e33SBarry Smith       if (sp) PetscCall(MatSetNullSpace(jac->schur, sp));
85320252d06SBarry Smith 
8549566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname));
8559566063dSJacob Faibussowitsch       PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, schurtestoption, NULL, &flg));
856514bf10dSMatthew G Knepley       if (flg) {
857514bf10dSMatthew G Knepley         DM  dmInner;
85821635b76SJed Brown         KSP kspInner;
859686bed4dSStefano Zampini         PC  pcInner;
86021635b76SJed Brown 
8619566063dSJacob Faibussowitsch         PetscCall(MatSchurComplementGetKSP(jac->schur, &kspInner));
8629566063dSJacob Faibussowitsch         PetscCall(KSPReset(kspInner));
8639566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(kspInner, jac->mat[0], jac->pmat[0]));
8649566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname));
86521635b76SJed Brown         /* Indent this deeper to emphasize the "inner" nature of this solver. */
8669566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject)pc, 2));
8679566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspInner->pc, (PetscObject)pc, 2));
8689566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(kspInner, schurprefix));
869514bf10dSMatthew G Knepley 
870514bf10dSMatthew G Knepley         /* Set DM for new solver */
8719566063dSJacob Faibussowitsch         PetscCall(KSPGetDM(jac->head->ksp, &dmInner));
8729566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(kspInner, dmInner));
8739566063dSJacob Faibussowitsch         PetscCall(KSPSetDMActive(kspInner, PETSC_FALSE));
874686bed4dSStefano Zampini 
875686bed4dSStefano Zampini         /* Defaults to PCKSP as preconditioner */
8769566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(kspInner, &pcInner));
8779566063dSJacob Faibussowitsch         PetscCall(PCSetType(pcInner, PCKSP));
8789566063dSJacob Faibussowitsch         PetscCall(PCKSPSetKSP(pcInner, jac->head->ksp));
879514bf10dSMatthew G Knepley       } else {
88021635b76SJed Brown         /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
88121635b76SJed Brown           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
88221635b76SJed Brown           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
88321635b76SJed 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
88421635b76SJed Brown           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
88521635b76SJed Brown           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
8869566063dSJacob Faibussowitsch         PetscCall(KSPSetType(jac->head->ksp, KSPGMRES));
8879566063dSJacob Faibussowitsch         PetscCall(MatSchurComplementSetKSP(jac->schur, jac->head->ksp));
888bafc1b83SMatthew G Knepley       }
8899566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->head->ksp, jac->mat[0], jac->pmat[0]));
8909566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(jac->head->ksp));
8919566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(jac->schur));
8923b224e63SBarry Smith 
8939566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)jac->schur, MATSCHURCOMPLEMENT, &flg));
894686bed4dSStefano Zampini       if (flg) { /* Need to do this otherwise PCSetUp_KSP will overwrite the amat of jac->head->ksp */
895686bed4dSStefano Zampini         KSP kspInner;
896686bed4dSStefano Zampini         PC  pcInner;
897686bed4dSStefano Zampini 
8989566063dSJacob Faibussowitsch         PetscCall(MatSchurComplementGetKSP(jac->schur, &kspInner));
8999566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(kspInner, &pcInner));
9009566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)pcInner, PCKSP, &flg));
901686bed4dSStefano Zampini         if (flg) {
902686bed4dSStefano Zampini           KSP ksp;
903686bed4dSStefano Zampini 
9049566063dSJacob Faibussowitsch           PetscCall(PCKSPGetKSP(pcInner, &ksp));
90548a46eb9SPierre Jolivet           if (ksp == jac->head->ksp) PetscCall(PCSetUseAmat(pcInner, PETSC_TRUE));
906686bed4dSStefano Zampini         }
907686bed4dSStefano Zampini       }
9089566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname));
9099566063dSJacob Faibussowitsch       PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, schurtestoption, NULL, &flg));
910443836d0SMatthew G Knepley       if (flg) {
911443836d0SMatthew G Knepley         DM dmInner;
912443836d0SMatthew G Knepley 
9139566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname));
9149566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper));
9159566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(jac->kspupper, pc->erroriffailure));
9169566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(jac->kspupper, schurprefix));
9179566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspupper, (PetscObject)pc, 1));
9189566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspupper->pc, (PetscObject)pc, 1));
9199566063dSJacob Faibussowitsch         PetscCall(KSPGetDM(jac->head->ksp, &dmInner));
9209566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(jac->kspupper, dmInner));
9219566063dSJacob Faibussowitsch         PetscCall(KSPSetDMActive(jac->kspupper, PETSC_FALSE));
9229566063dSJacob Faibussowitsch         PetscCall(KSPSetFromOptions(jac->kspupper));
9239566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(jac->kspupper, jac->mat[0], jac->pmat[0]));
9249566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(jac->head->x, &jac->head->z));
925443836d0SMatthew G Knepley       } else {
926443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
9279566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)jac->head->ksp));
928443836d0SMatthew G Knepley       }
929443836d0SMatthew G Knepley 
93048a46eb9SPierre Jolivet       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) PetscCall(MatSchurComplementGetPmat(jac->schur, MAT_INITIAL_MATRIX, &jac->schurp));
9319566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspschur));
9329566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(jac->kspschur, pc->erroriffailure));
9339566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspschur, (PetscObject)pc, 1));
934084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
9357233a360SDmitry Karpeev         PC pcschur;
9369566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(jac->kspschur, &pcschur));
9379566063dSJacob Faibussowitsch         PetscCall(PCSetType(pcschur, PCNONE));
938084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
939e74569cdSMatthew G. Knepley       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
9409566063dSJacob Faibussowitsch         PetscCall(MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user));
941e69d4d44SBarry Smith       }
9429566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->kspschur, jac->schur, FieldSplitSchurPre(jac)));
9439566063dSJacob Faibussowitsch       PetscCall(KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix));
9449566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(jac->kspschur, Dprefix));
945c096484dSStefano Zampini       /* propagate DM */
946b20b4189SMatthew G. Knepley       {
947b20b4189SMatthew G. Knepley         DM sdm;
9489566063dSJacob Faibussowitsch         PetscCall(KSPGetDM(jac->head->next->ksp, &sdm));
949b20b4189SMatthew G. Knepley         if (sdm) {
9509566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(jac->kspschur, sdm));
9519566063dSJacob Faibussowitsch           PetscCall(KSPSetDMActive(jac->kspschur, PETSC_FALSE));
952b20b4189SMatthew G. Knepley         }
953b20b4189SMatthew G. Knepley       }
9543b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
95520b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
9569566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(jac->kspschur));
9573b224e63SBarry Smith     }
9589566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac->schur, MAT_FINAL_ASSEMBLY));
9599566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac->schur, MAT_FINAL_ASSEMBLY));
960093c86ffSJed Brown 
961093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
9629566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(lscname, sizeof(lscname), "%s_LSC_L", ilink->splitname));
9639566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)pc->mat, lscname, (PetscObject *)&LSC_L));
9649566063dSJacob Faibussowitsch     if (!LSC_L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, lscname, (PetscObject *)&LSC_L));
9659566063dSJacob Faibussowitsch     if (LSC_L) PetscCall(PetscObjectCompose((PetscObject)jac->schur, "LSC_L", (PetscObject)LSC_L));
9669566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(lscname, sizeof(lscname), "%s_LSC_Lp", ilink->splitname));
9679566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)pc->pmat, lscname, (PetscObject *)&LSC_L));
9689566063dSJacob Faibussowitsch     if (!LSC_L) PetscCall(PetscObjectQuery((PetscObject)pc->mat, lscname, (PetscObject *)&LSC_L));
9699566063dSJacob Faibussowitsch     if (LSC_L) PetscCall(PetscObjectCompose((PetscObject)jac->schur, "LSC_Lp", (PetscObject)LSC_L));
970a51937d4SCarola Kruse   } else if (jac->type == PC_COMPOSITE_GKB) {
971a51937d4SCarola Kruse     IS       ccis;
972a51937d4SCarola Kruse     PetscInt rstart, rend;
973a51937d4SCarola Kruse 
97408401ef6SPierre Jolivet     PetscCheck(nsplit == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "To use GKB preconditioner you must have exactly 2 fields");
975a51937d4SCarola Kruse 
976a51937d4SCarola Kruse     ilink = jac->head;
977a51937d4SCarola Kruse 
978a51937d4SCarola Kruse     /* When extracting off-diagonal submatrices, we take complements from this range */
9799566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(pc->mat, &rstart, &rend));
980a51937d4SCarola Kruse 
9819566063dSJacob Faibussowitsch     PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis));
982a51937d4SCarola Kruse     if (jac->offdiag_use_amat) {
9839566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B));
984a51937d4SCarola Kruse     } else {
9859566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B));
986a51937d4SCarola Kruse     }
9879566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ccis));
988e071a0a4SCarola Kruse     /* Create work vectors for GKB algorithm */
9899566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ilink->x, &jac->u));
9909566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ilink->x, &jac->Hu));
9919566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ilink->x, &jac->w2));
992a51937d4SCarola Kruse     ilink = ilink->next;
9939566063dSJacob Faibussowitsch     PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis));
994a51937d4SCarola Kruse     if (jac->offdiag_use_amat) {
9959566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C));
996a51937d4SCarola Kruse     } else {
9979566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C));
998a51937d4SCarola Kruse     }
9999566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ccis));
1000e071a0a4SCarola Kruse     /* Create work vectors for GKB algorithm */
10019566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ilink->x, &jac->v));
10029566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ilink->x, &jac->d));
10039566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ilink->x, &jac->w1));
10049566063dSJacob Faibussowitsch     PetscCall(MatGolubKahanComputeExplicitOperator(jac->mat[0], jac->B, jac->C, &jac->H, jac->gkbnu));
10059566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(jac->gkbdelay, &jac->vecz));
1006e071a0a4SCarola Kruse 
1007a51937d4SCarola Kruse     ilink = jac->head;
10089566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ilink->ksp, jac->H, jac->H));
10099566063dSJacob Faibussowitsch     if (!jac->suboptionsset) PetscCall(KSPSetFromOptions(ilink->ksp));
1010e071a0a4SCarola Kruse     /* Create gkb_monitor context */
1011de482cd7SCarola Kruse     if (jac->gkbmonitor) {
1012de482cd7SCarola Kruse       PetscInt tablevel;
10139566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &jac->gkbviewer));
10149566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(jac->gkbviewer, PETSCVIEWERASCII));
10159566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetTabLevel((PetscObject)ilink->ksp, &tablevel));
10169566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISetTab(jac->gkbviewer, tablevel));
10179566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->ksp, (PetscObject)ilink->ksp, 1));
1018de482cd7SCarola Kruse     }
10193b224e63SBarry Smith   } else {
102068bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
102197bbdb24SBarry Smith     i     = 0;
10225a9f2f41SSatish Balay     ilink = jac->head;
10235a9f2f41SSatish Balay     while (ilink) {
10249566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(ilink->ksp, jac->mat[i], jac->pmat[i]));
10253b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
10269566063dSJacob Faibussowitsch       if (!jac->suboptionsset) PetscCall(KSPSetFromOptions(ilink->ksp));
102797bbdb24SBarry Smith       i++;
10285a9f2f41SSatish Balay       ilink = ilink->next;
10290971522cSBarry Smith     }
10303b224e63SBarry Smith   }
10313b224e63SBarry Smith 
10325ddf11f8SNicolas Barnafi   /* Set coordinates to the sub PC objects whenever these are set */
10335ddf11f8SNicolas Barnafi   if (jac->coordinates_set) {
10345ddf11f8SNicolas Barnafi     PC pc_coords;
10355ddf11f8SNicolas Barnafi     if (jac->type == PC_COMPOSITE_SCHUR) {
10365ddf11f8SNicolas Barnafi       // Head is first block.
10379566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(jac->head->ksp, &pc_coords));
10389566063dSJacob Faibussowitsch       PetscCall(PCSetCoordinates(pc_coords, jac->head->dim, jac->head->ndofs, jac->head->coords));
10395ddf11f8SNicolas Barnafi       // Second one is Schur block, but its KSP object is in kspschur.
10409566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(jac->kspschur, &pc_coords));
10419566063dSJacob Faibussowitsch       PetscCall(PCSetCoordinates(pc_coords, jac->head->next->dim, jac->head->next->ndofs, jac->head->next->coords));
10425ddf11f8SNicolas Barnafi     } else if (jac->type == PC_COMPOSITE_GKB) {
10439566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Warning: Setting coordinates does nothing for the GKB Fieldpslit preconditioner"));
10445ddf11f8SNicolas Barnafi     } else {
10455ddf11f8SNicolas Barnafi       ilink = jac->head;
10465ddf11f8SNicolas Barnafi       while (ilink) {
10479566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ilink->ksp, &pc_coords));
10489566063dSJacob Faibussowitsch         PetscCall(PCSetCoordinates(pc_coords, ilink->dim, ilink->ndofs, ilink->coords));
10495ddf11f8SNicolas Barnafi         ilink = ilink->next;
10505ddf11f8SNicolas Barnafi       }
10515ddf11f8SNicolas Barnafi     }
10525ddf11f8SNicolas Barnafi   }
10535ddf11f8SNicolas Barnafi 
1054c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
10550971522cSBarry Smith   PetscFunctionReturn(0);
10560971522cSBarry Smith }
10570971522cSBarry Smith 
10585a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink, xx, yy) \
10599371c9d4SSatish Balay   (VecScatterBegin(ilink->sctx, xx, ilink->x, INSERT_VALUES, SCATTER_FORWARD) || VecScatterEnd(ilink->sctx, xx, ilink->x, INSERT_VALUES, SCATTER_FORWARD) || PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL) || \
10609371c9d4SSatish Balay    KSPSolve(ilink->ksp, ilink->x, ilink->y) || KSPCheckSolve(ilink->ksp, pc, ilink->y) || PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL) || VecScatterBegin(ilink->sctx, ilink->y, yy, ADD_VALUES, SCATTER_REVERSE) || \
1061ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx, ilink->y, yy, ADD_VALUES, SCATTER_REVERSE))
106279416396SBarry Smith 
10639371c9d4SSatish Balay static PetscErrorCode PCApply_FieldSplit_Schur(PC pc, Vec x, Vec y) {
10643b224e63SBarry Smith   PC_FieldSplit    *jac    = (PC_FieldSplit *)pc->data;
10653b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
1066443836d0SMatthew G Knepley   KSP               kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
10673b224e63SBarry Smith 
10683b224e63SBarry Smith   PetscFunctionBegin;
1069c5d2311dSJed Brown   switch (jac->schurfactorization) {
1070c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
1071a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
10729566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
10739566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD));
10749566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
10759566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
10769566063dSJacob Faibussowitsch     PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y));
10779566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y));
10789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
10799566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
10809566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD));
10819566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
10829566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y));
10839566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y));
10849566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
10859566063dSJacob Faibussowitsch     PetscCall(VecScale(ilinkD->y, jac->schurscale));
10869566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
10879566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
10889566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
1089c5d2311dSJed Brown     break;
1090c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
1091a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
10929566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
10939566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
10949566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
10959566063dSJacob Faibussowitsch     PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y));
10969566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y));
10979566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
10989566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->C, ilinkA->y, ilinkD->x));
10999566063dSJacob Faibussowitsch     PetscCall(VecScale(ilinkD->x, -1.));
11009566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD));
11019566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
11029566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD));
11039566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
11049566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y));
11059566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y));
11069566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
11079566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
11089566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
11099566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
1110c5d2311dSJed Brown     break;
1111c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
1112a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
11139566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD));
11149566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD));
11159566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
11169566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y));
11179566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y));
11189371c9d4SSatish Balay     PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
11199371c9d4SSatish Balay     PetscCall(MatMult(jac->B, ilinkD->y, ilinkA->x));
11209566063dSJacob Faibussowitsch     PetscCall(VecScale(ilinkA->x, -1.));
11219566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, ADD_VALUES, SCATTER_FORWARD));
11229566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
11239566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, ADD_VALUES, SCATTER_FORWARD));
11249566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
11259566063dSJacob Faibussowitsch     PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y));
11269566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y));
11279566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
11289566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
11299566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
11309566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
1131c5d2311dSJed Brown     break;
1132c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
1133c238f8cdSStefano Zampini     /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1] */
11349566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
11359566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
11369566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(KSP_Solve_FS_L, kspLower, ilinkA->x, ilinkA->y, NULL));
11379566063dSJacob Faibussowitsch     PetscCall(KSPSolve(kspLower, ilinkA->x, ilinkA->y));
11389566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(kspLower, pc, ilinkA->y));
11399566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(KSP_Solve_FS_L, kspLower, ilinkA->x, ilinkA->y, NULL));
11409566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->C, ilinkA->y, ilinkD->x));
11419566063dSJacob Faibussowitsch     PetscCall(VecScale(ilinkD->x, -1.0));
11429566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD));
11439566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD));
11443b224e63SBarry Smith 
11459566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
11469566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y));
11479566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y));
11489566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL));
11499566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
11503b224e63SBarry Smith 
1151443836d0SMatthew G Knepley     if (kspUpper == kspA) {
11529566063dSJacob Faibussowitsch       PetscCall(MatMult(jac->B, ilinkD->y, ilinkA->y));
11539566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ilinkA->x, -1.0, ilinkA->y));
11549566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
11559566063dSJacob Faibussowitsch       PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y));
11569566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y));
11579566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
1158443836d0SMatthew G Knepley     } else {
11599566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL));
11609566063dSJacob Faibussowitsch       PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y));
11619566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y));
11629566063dSJacob Faibussowitsch       PetscCall(MatMult(jac->B, ilinkD->y, ilinkA->x));
11639566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(KSP_Solve_FS_U, kspUpper, ilinkA->x, ilinkA->z, NULL));
11649566063dSJacob Faibussowitsch       PetscCall(KSPSolve(kspUpper, ilinkA->x, ilinkA->z));
11659566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(kspUpper, pc, ilinkA->z));
11669566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(KSP_Solve_FS_U, kspUpper, ilinkA->x, ilinkA->z, NULL));
11679566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ilinkA->y, -1.0, ilinkA->z));
1168443836d0SMatthew G Knepley     }
11699566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
11709566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
11719566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
1172c5d2311dSJed Brown   }
11733b224e63SBarry Smith   PetscFunctionReturn(0);
11743b224e63SBarry Smith }
11753b224e63SBarry Smith 
11769371c9d4SSatish Balay static PetscErrorCode PCApply_FieldSplit(PC pc, Vec x, Vec y) {
11770971522cSBarry Smith   PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
11785a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
1179939b8a20SBarry Smith   PetscInt          cnt, bs;
11800971522cSBarry Smith 
11810971522cSBarry Smith   PetscFunctionBegin;
118279416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
11831b9fc7fcSBarry Smith     if (jac->defaultsplit) {
11849566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(x, &bs));
11852472a847SBarry Smith       PetscCheck(jac->bs <= 0 || bs == jac->bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Blocksize of x vector %" PetscInt_FMT " does not match fieldsplit blocksize %" PetscInt_FMT, bs, jac->bs);
11869566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(y, &bs));
11872472a847SBarry Smith       PetscCheck(jac->bs <= 0 || bs == jac->bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Blocksize of y vector %" PetscInt_FMT " does not match fieldsplit blocksize %" PetscInt_FMT, bs, jac->bs);
11889566063dSJacob Faibussowitsch       PetscCall(VecStrideGatherAll(x, jac->x, INSERT_VALUES));
11895a9f2f41SSatish Balay       while (ilink) {
11909566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
11919566063dSJacob Faibussowitsch         PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y));
11929566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y));
11939566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
11945a9f2f41SSatish Balay         ilink = ilink->next;
11950971522cSBarry Smith       }
11969566063dSJacob Faibussowitsch       PetscCall(VecStrideScatterAll(jac->y, y, INSERT_VALUES));
11971b9fc7fcSBarry Smith     } else {
11989566063dSJacob Faibussowitsch       PetscCall(VecSet(y, 0.0));
11995a9f2f41SSatish Balay       while (ilink) {
12009566063dSJacob Faibussowitsch         PetscCall(FieldSplitSplitSolveAdd(ilink, x, y));
12015a9f2f41SSatish Balay         ilink = ilink->next;
12021b9fc7fcSBarry Smith       }
12031b9fc7fcSBarry Smith     }
1204e52d2c62SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
12059566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
1206e52d2c62SBarry Smith     /* solve on first block for first block variables */
12079566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, INSERT_VALUES, SCATTER_FORWARD));
12089566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, INSERT_VALUES, SCATTER_FORWARD));
12099566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12109566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y));
12119566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y));
12129566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12139566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
12149566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
1215e52d2c62SBarry Smith 
1216e52d2c62SBarry Smith     /* compute the residual only onto second block variables using first block variables */
12179566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->Afield[1], ilink->y, ilink->next->x));
1218e52d2c62SBarry Smith     ilink = ilink->next;
12199566063dSJacob Faibussowitsch     PetscCall(VecScale(ilink->x, -1.0));
12209566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD));
12219566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD));
1222e52d2c62SBarry Smith 
1223e52d2c62SBarry Smith     /* solve on second block variables */
12249566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12259566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y));
12269566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y));
12279566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12289566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
12299566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
123016913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
123179416396SBarry Smith     if (!jac->w1) {
12329566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &jac->w1));
12339566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &jac->w2));
123479416396SBarry Smith     }
12359566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
12369566063dSJacob Faibussowitsch     PetscCall(FieldSplitSplitSolveAdd(ilink, x, y));
12373e197d65SBarry Smith     cnt = 1;
12385a9f2f41SSatish Balay     while (ilink->next) {
12395a9f2f41SSatish Balay       ilink = ilink->next;
12403e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
12419566063dSJacob Faibussowitsch       PetscCall(MatMult(jac->Afield[cnt++], y, ilink->x));
12429566063dSJacob Faibussowitsch       PetscCall(VecScale(ilink->x, -1.0));
12439566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD));
12449566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD));
12459566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12469566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y));
12479566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y));
12489566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12499566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
12509566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
12513e197d65SBarry Smith     }
125251f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
125311755939SBarry Smith       cnt -= 2;
125451f519a2SBarry Smith       while (ilink->previous) {
125551f519a2SBarry Smith         ilink = ilink->previous;
125611755939SBarry Smith         /* compute the residual only over the part of the vector needed */
12579566063dSJacob Faibussowitsch         PetscCall(MatMult(jac->Afield[cnt--], y, ilink->x));
12589566063dSJacob Faibussowitsch         PetscCall(VecScale(ilink->x, -1.0));
12599566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD));
12609566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD));
12619566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12629566063dSJacob Faibussowitsch         PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y));
12639566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y));
12649566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
12659566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
12669566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE));
126751f519a2SBarry Smith       }
126811755939SBarry Smith     }
126963a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unsupported or unknown composition %d", (int)jac->type);
12700971522cSBarry Smith   PetscFunctionReturn(0);
12710971522cSBarry Smith }
12720971522cSBarry Smith 
12739371c9d4SSatish Balay static PetscErrorCode PCApply_FieldSplit_GKB(PC pc, Vec x, Vec y) {
1274a51937d4SCarola Kruse   PC_FieldSplit    *jac    = (PC_FieldSplit *)pc->data;
1275a51937d4SCarola Kruse   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
1276a51937d4SCarola Kruse   KSP               ksp = ilinkA->ksp;
1277de482cd7SCarola Kruse   Vec               u, v, Hu, d, work1, work2;
1278e071a0a4SCarola Kruse   PetscScalar       alpha, z, nrmz2, *vecz;
1279e071a0a4SCarola Kruse   PetscReal         lowbnd, nu, beta;
1280a51937d4SCarola Kruse   PetscInt          j, iterGKB;
1281a51937d4SCarola Kruse 
1282a51937d4SCarola Kruse   PetscFunctionBegin;
12839566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
12849566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD));
12859566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD));
12869566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD));
1287e071a0a4SCarola Kruse 
1288e071a0a4SCarola Kruse   u     = jac->u;
1289e071a0a4SCarola Kruse   v     = jac->v;
1290e071a0a4SCarola Kruse   Hu    = jac->Hu;
1291e071a0a4SCarola Kruse   d     = jac->d;
1292e071a0a4SCarola Kruse   work1 = jac->w1;
1293e071a0a4SCarola Kruse   work2 = jac->w2;
1294e071a0a4SCarola Kruse   vecz  = jac->vecz;
1295a51937d4SCarola Kruse 
1296a51937d4SCarola Kruse   /* Change RHS to comply with matrix regularization H = A + nu*B*B' */
1297a51937d4SCarola Kruse   /* Add q = q + nu*B*b */
1298a51937d4SCarola Kruse   if (jac->gkbnu) {
1299a51937d4SCarola Kruse     nu = jac->gkbnu;
13009566063dSJacob Faibussowitsch     PetscCall(VecScale(ilinkD->x, jac->gkbnu));
13019566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(jac->B, ilinkD->x, ilinkA->x, ilinkA->x)); /* q = q + nu*B*b */
1302a51937d4SCarola Kruse   } else {
1303a51937d4SCarola Kruse     /* Situation when no augmented Lagrangian is used. Then we set inner  */
1304a51937d4SCarola Kruse     /* matrix N = I in [Ar13], and thus nu = 1.                           */
1305a51937d4SCarola Kruse     nu = 1;
1306a51937d4SCarola Kruse   }
1307a51937d4SCarola Kruse 
1308a51937d4SCarola Kruse   /* Transform rhs from [q,tilde{b}] to [0,b] */
13099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(ilinkA->event, ksp, ilinkA->x, ilinkA->y, NULL));
13109566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, ilinkA->x, ilinkA->y));
13119566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(ksp, pc, ilinkA->y));
13129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(ilinkA->event, ksp, ilinkA->x, ilinkA->y, NULL));
13139566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(jac->B, ilinkA->y, work1));
13149566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(work1, 1.0 / nu, -1.0, ilinkD->x)); /* c = b - B'*x        */
1315a51937d4SCarola Kruse 
1316a51937d4SCarola Kruse   /* First step of algorithm */
13179566063dSJacob Faibussowitsch   PetscCall(VecNorm(work1, NORM_2, &beta)); /* beta = sqrt(nu*c'*c)*/
1318e071a0a4SCarola Kruse   KSPCheckDot(ksp, beta);
1319addd1e01SJunchao Zhang   beta = PetscSqrtReal(nu) * beta;
13209566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(v, nu / beta, 0.0, work1)); /* v = nu/beta *c      */
13219566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->B, v, work2));          /* u = H^{-1}*B*v      */
13229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(ilinkA->event, ksp, work2, u, NULL));
13239566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, work2, u));
13249566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(ksp, pc, u));
13259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(ilinkA->event, ksp, work2, u, NULL));
13269566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->H, u, Hu)); /* alpha = u'*H*u      */
13279566063dSJacob Faibussowitsch   PetscCall(VecDot(Hu, u, &alpha));
1328e071a0a4SCarola Kruse   KSPCheckDot(ksp, alpha);
132908401ef6SPierre Jolivet   PetscCheck(PetscRealPart(alpha) > 0.0, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "GKB preconditioner diverged, H is not positive definite");
1330addd1e01SJunchao Zhang   alpha = PetscSqrtReal(PetscAbsScalar(alpha));
13319566063dSJacob Faibussowitsch   PetscCall(VecScale(u, 1.0 / alpha));
13329566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(d, 1.0 / alpha, 0.0, v)); /* v = nu/beta *c      */
1333de482cd7SCarola Kruse 
1334a51937d4SCarola Kruse   z       = beta / alpha;
1335a51937d4SCarola Kruse   vecz[1] = z;
1336a51937d4SCarola Kruse 
1337de482cd7SCarola Kruse   /* Computation of first iterate x(1) and p(1) */
13389566063dSJacob Faibussowitsch   PetscCall(VecAXPY(ilinkA->y, z, u));
13399566063dSJacob Faibussowitsch   PetscCall(VecCopy(d, ilinkD->y));
13409566063dSJacob Faibussowitsch   PetscCall(VecScale(ilinkD->y, -z));
1341a51937d4SCarola Kruse 
13429371c9d4SSatish Balay   iterGKB = 1;
13439371c9d4SSatish Balay   lowbnd  = 2 * jac->gkbtol;
134448a46eb9SPierre Jolivet   if (jac->gkbmonitor) PetscCall(PetscViewerASCIIPrintf(jac->gkbviewer, "%3" PetscInt_FMT " GKB Lower bound estimate %14.12e\n", iterGKB, (double)lowbnd));
1345de482cd7SCarola Kruse 
1346a51937d4SCarola Kruse   while (iterGKB < jac->gkbmaxit && lowbnd > jac->gkbtol) {
1347a51937d4SCarola Kruse     iterGKB += 1;
13489566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(jac->B, u, work1)); /* v <- nu*(B'*u-alpha/nu*v) */
13499566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(v, nu, -alpha, work1));
13509566063dSJacob Faibussowitsch     PetscCall(VecNorm(v, NORM_2, &beta)); /* beta = sqrt(nu)*v'*v      */
1351addd1e01SJunchao Zhang     beta = beta / PetscSqrtReal(nu);
13529566063dSJacob Faibussowitsch     PetscCall(VecScale(v, 1.0 / beta));
13539566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->B, v, work2)); /* u <- H^{-1}*(B*v-beta*H*u) */
13549566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->H, u, Hu));
13559566063dSJacob Faibussowitsch     PetscCall(VecAXPY(work2, -beta, Hu));
13569566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ilinkA->event, ksp, work2, u, NULL));
13579566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ksp, work2, u));
13589566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(ksp, pc, u));
13599566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ilinkA->event, ksp, work2, u, NULL));
13609566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->H, u, Hu)); /* alpha = u'*H*u            */
13619566063dSJacob Faibussowitsch     PetscCall(VecDot(Hu, u, &alpha));
1362e071a0a4SCarola Kruse     KSPCheckDot(ksp, alpha);
136308401ef6SPierre Jolivet     PetscCheck(PetscRealPart(alpha) > 0.0, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "GKB preconditioner diverged, H is not positive definite");
1364addd1e01SJunchao Zhang     alpha = PetscSqrtReal(PetscAbsScalar(alpha));
13659566063dSJacob Faibussowitsch     PetscCall(VecScale(u, 1.0 / alpha));
1366a51937d4SCarola Kruse 
1367e071a0a4SCarola Kruse     z       = -beta / alpha * z; /* z <- beta/alpha*z     */
1368a51937d4SCarola Kruse     vecz[0] = z;
1369a51937d4SCarola Kruse 
1370a51937d4SCarola Kruse     /* Computation of new iterate x(i+1) and p(i+1) */
13719566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(d, 1.0 / alpha, -beta / alpha, v)); /* d = (v-beta*d)/alpha */
13729566063dSJacob Faibussowitsch     PetscCall(VecAXPY(ilinkA->y, z, u));                   /* r = r + z*u          */
13739566063dSJacob Faibussowitsch     PetscCall(VecAXPY(ilinkD->y, -z, d));                  /* p = p - z*d          */
13749566063dSJacob Faibussowitsch     PetscCall(MatMult(jac->H, ilinkA->y, Hu));             /* ||u||_H = u'*H*u     */
13759566063dSJacob Faibussowitsch     PetscCall(VecDot(Hu, ilinkA->y, &nrmz2));
1376a51937d4SCarola Kruse 
1377a51937d4SCarola Kruse     /* Compute Lower Bound estimate */
1378a51937d4SCarola Kruse     if (iterGKB > jac->gkbdelay) {
1379a51937d4SCarola Kruse       lowbnd = 0.0;
1380ad540459SPierre Jolivet       for (j = 0; j < jac->gkbdelay; j++) lowbnd += PetscAbsScalar(vecz[j] * vecz[j]);
1381addd1e01SJunchao Zhang       lowbnd = PetscSqrtReal(lowbnd / PetscAbsScalar(nrmz2));
1382a51937d4SCarola Kruse     }
1383a51937d4SCarola Kruse 
1384ad540459SPierre Jolivet     for (j = 0; j < jac->gkbdelay - 1; j++) vecz[jac->gkbdelay - j - 1] = vecz[jac->gkbdelay - j - 2];
138548a46eb9SPierre Jolivet     if (jac->gkbmonitor) PetscCall(PetscViewerASCIIPrintf(jac->gkbviewer, "%3" PetscInt_FMT " GKB Lower bound estimate %14.12e\n", iterGKB, (double)lowbnd));
1386a51937d4SCarola Kruse   }
1387a51937d4SCarola Kruse 
13889566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
13899566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE));
13909566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
13919566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE));
1392a51937d4SCarola Kruse 
1393a51937d4SCarola Kruse   PetscFunctionReturn(0);
1394a51937d4SCarola Kruse }
1395a51937d4SCarola Kruse 
1396421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink, xx, yy) \
13979371c9d4SSatish Balay   (VecScatterBegin(ilink->sctx, xx, ilink->y, INSERT_VALUES, SCATTER_FORWARD) || VecScatterEnd(ilink->sctx, xx, ilink->y, INSERT_VALUES, SCATTER_FORWARD) || PetscLogEventBegin(ilink->event, ilink->ksp, ilink->y, ilink->x, NULL) || \
13989371c9d4SSatish Balay    KSPSolveTranspose(ilink->ksp, ilink->y, ilink->x) || KSPCheckSolve(ilink->ksp, pc, ilink->x) || PetscLogEventEnd(ilink->event, ilink->ksp, ilink->y, ilink->x, NULL) || VecScatterBegin(ilink->sctx, ilink->x, yy, ADD_VALUES, SCATTER_REVERSE) || \
1399ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx, ilink->x, yy, ADD_VALUES, SCATTER_REVERSE))
1400421e10b8SBarry Smith 
14019371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc, Vec x, Vec y) {
1402421e10b8SBarry Smith   PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
1403421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1404939b8a20SBarry Smith   PetscInt          bs;
1405421e10b8SBarry Smith 
1406421e10b8SBarry Smith   PetscFunctionBegin;
1407421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1408421e10b8SBarry Smith     if (jac->defaultsplit) {
14099566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(x, &bs));
14102472a847SBarry Smith       PetscCheck(jac->bs <= 0 || bs == jac->bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Blocksize of x vector %" PetscInt_FMT " does not match fieldsplit blocksize %" PetscInt_FMT, bs, jac->bs);
14119566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(y, &bs));
14122472a847SBarry Smith       PetscCheck(jac->bs <= 0 || bs == jac->bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Blocksize of y vector %" PetscInt_FMT " does not match fieldsplit blocksize %" PetscInt_FMT, bs, jac->bs);
14139566063dSJacob Faibussowitsch       PetscCall(VecStrideGatherAll(x, jac->x, INSERT_VALUES));
1414421e10b8SBarry Smith       while (ilink) {
14159566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
14169566063dSJacob Faibussowitsch         PetscCall(KSPSolveTranspose(ilink->ksp, ilink->x, ilink->y));
14179566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y));
14189566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL));
1419421e10b8SBarry Smith         ilink = ilink->next;
1420421e10b8SBarry Smith       }
14219566063dSJacob Faibussowitsch       PetscCall(VecStrideScatterAll(jac->y, y, INSERT_VALUES));
1422421e10b8SBarry Smith     } else {
14239566063dSJacob Faibussowitsch       PetscCall(VecSet(y, 0.0));
1424421e10b8SBarry Smith       while (ilink) {
14259566063dSJacob Faibussowitsch         PetscCall(FieldSplitSplitSolveAddTranspose(ilink, x, y));
1426421e10b8SBarry Smith         ilink = ilink->next;
1427421e10b8SBarry Smith       }
1428421e10b8SBarry Smith     }
1429421e10b8SBarry Smith   } else {
1430421e10b8SBarry Smith     if (!jac->w1) {
14319566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &jac->w1));
14329566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &jac->w2));
1433421e10b8SBarry Smith     }
14349566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
1435421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
14369566063dSJacob Faibussowitsch       PetscCall(FieldSplitSplitSolveAddTranspose(ilink, x, y));
1437421e10b8SBarry Smith       while (ilink->next) {
1438421e10b8SBarry Smith         ilink = ilink->next;
14399566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(pc->mat, y, jac->w1));
14409566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(jac->w2, -1.0, jac->w1, x));
14419566063dSJacob Faibussowitsch         PetscCall(FieldSplitSplitSolveAddTranspose(ilink, jac->w2, y));
1442421e10b8SBarry Smith       }
1443421e10b8SBarry Smith       while (ilink->previous) {
1444421e10b8SBarry Smith         ilink = ilink->previous;
14459566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(pc->mat, y, jac->w1));
14469566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(jac->w2, -1.0, jac->w1, x));
14479566063dSJacob Faibussowitsch         PetscCall(FieldSplitSplitSolveAddTranspose(ilink, jac->w2, y));
1448421e10b8SBarry Smith       }
1449421e10b8SBarry Smith     } else {
1450421e10b8SBarry Smith       while (ilink->next) { /* get to last entry in linked list */
1451421e10b8SBarry Smith         ilink = ilink->next;
1452421e10b8SBarry Smith       }
14539566063dSJacob Faibussowitsch       PetscCall(FieldSplitSplitSolveAddTranspose(ilink, x, y));
1454421e10b8SBarry Smith       while (ilink->previous) {
1455421e10b8SBarry Smith         ilink = ilink->previous;
14569566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(pc->mat, y, jac->w1));
14579566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(jac->w2, -1.0, jac->w1, x));
14589566063dSJacob Faibussowitsch         PetscCall(FieldSplitSplitSolveAddTranspose(ilink, jac->w2, y));
1459421e10b8SBarry Smith       }
1460421e10b8SBarry Smith     }
1461421e10b8SBarry Smith   }
1462421e10b8SBarry Smith   PetscFunctionReturn(0);
1463421e10b8SBarry Smith }
1464421e10b8SBarry Smith 
14659371c9d4SSatish Balay static PetscErrorCode PCReset_FieldSplit(PC pc) {
14660971522cSBarry Smith   PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
14675a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head, next;
14680971522cSBarry Smith 
14690971522cSBarry Smith   PetscFunctionBegin;
14705a9f2f41SSatish Balay   while (ilink) {
14719566063dSJacob Faibussowitsch     PetscCall(KSPDestroy(&ilink->ksp));
14729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ilink->x));
14739566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ilink->y));
14749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ilink->z));
14759566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&ilink->sctx));
14769566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ilink->is));
14779566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ilink->is_col));
14789566063dSJacob Faibussowitsch     PetscCall(PetscFree(ilink->splitname));
14799566063dSJacob Faibussowitsch     PetscCall(PetscFree(ilink->fields));
14809566063dSJacob Faibussowitsch     PetscCall(PetscFree(ilink->fields_col));
14815a9f2f41SSatish Balay     next = ilink->next;
14829566063dSJacob Faibussowitsch     PetscCall(PetscFree(ilink));
14835a9f2f41SSatish Balay     ilink = next;
14840971522cSBarry Smith   }
1485f5f0d762SBarry Smith   jac->head = NULL;
14869566063dSJacob Faibussowitsch   PetscCall(PetscFree2(jac->x, jac->y));
1487574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
14889566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->nsplits, &jac->mat));
1489574deadeSBarry Smith   } else if (jac->mat) {
14900298fd71SBarry Smith     jac->mat = NULL;
1491574deadeSBarry Smith   }
14929566063dSJacob Faibussowitsch   if (jac->pmat) PetscCall(MatDestroyMatrices(jac->nsplits, &jac->pmat));
14939566063dSJacob Faibussowitsch   if (jac->Afield) PetscCall(MatDestroyMatrices(jac->nsplits, &jac->Afield));
1494f5f0d762SBarry Smith   jac->nsplits = 0;
14959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->w1));
14969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->w2));
14979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->schur));
14989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->schurp));
14999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->schur_user));
15009566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->kspschur));
15019566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->kspupper));
15029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->B));
15039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->C));
15049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->H));
15059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->u));
15069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->v));
15079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->Hu));
15089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->d));
15099566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->vecz));
15109566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&jac->gkbviewer));
15116dbb499eSCian Wilson   jac->isrestrict = PETSC_FALSE;
1512574deadeSBarry Smith   PetscFunctionReturn(0);
1513574deadeSBarry Smith }
1514574deadeSBarry Smith 
15159371c9d4SSatish Balay static PetscErrorCode PCDestroy_FieldSplit(PC pc) {
1516574deadeSBarry Smith   PetscFunctionBegin;
15179566063dSJacob Faibussowitsch   PetscCall(PCReset_FieldSplit(pc));
15189566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
15192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
15209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetFields_C", NULL));
15219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetIS_C", NULL));
15229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetType_C", NULL));
15239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetBlockSize_C", NULL));
15242e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitRestrictIS_C", NULL));
15252e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSchurGetSubKSP_C", NULL));
15262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", NULL));
15272e956fe4SStefano Zampini 
15282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBTol_C", NULL));
15292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBMaxit_C", NULL));
15302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBNu_C", NULL));
15312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBDelay_C", NULL));
15322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", NULL));
15339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurPre_C", NULL));
15349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSchurPre_C", NULL));
15359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurFactType_C", NULL));
15362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurScale_C", NULL));
15370971522cSBarry Smith   PetscFunctionReturn(0);
15380971522cSBarry Smith }
15390971522cSBarry Smith 
15409371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc, PetscOptionItems *PetscOptionsObject) {
15416c924f48SJed Brown   PetscInt        bs;
15427b752e3dSPatrick Sanan   PetscBool       flg;
15439dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit *)pc->data;
15443b224e63SBarry Smith   PCCompositeType ctype;
15451b9fc7fcSBarry Smith 
15460971522cSBarry Smith   PetscFunctionBegin;
1547d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "FieldSplit options");
15489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_fieldsplit_dm_splits", "Whether to use DMCreateFieldDecomposition() for splits", "PCFieldSplitSetDMSplits", jac->dm_splits, &jac->dm_splits, NULL));
15499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_fieldsplit_block_size", "Blocksize that defines number of fields", "PCFieldSplitSetBlockSize", jac->bs, &bs, &flg));
15501baa6e33SBarry Smith   if (flg) PetscCall(PCFieldSplitSetBlockSize(pc, bs));
15512686e3e9SMatthew G. Knepley   jac->diag_use_amat = pc->useAmat;
15529566063dSJacob Faibussowitsch   PetscCall(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));
15532686e3e9SMatthew G. Knepley   jac->offdiag_use_amat = pc->useAmat;
15549566063dSJacob Faibussowitsch   PetscCall(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));
15559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_fieldsplit_detect_saddle_point", "Form 2-way split by detecting zero diagonal entries", "PCFieldSplitSetDetectSaddlePoint", jac->detect, &jac->detect, NULL));
15569566063dSJacob Faibussowitsch   PetscCall(PCFieldSplitSetDetectSaddlePoint(pc, jac->detect)); /* Sets split type and Schur PC type */
15579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_fieldsplit_type", "Type of composition", "PCFieldSplitSetType", PCCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&ctype, &flg));
15581baa6e33SBarry Smith   if (flg) PetscCall(PCFieldSplitSetType(pc, ctype));
1559c30613efSMatthew Knepley   /* Only setup fields once */
1560c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1561d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1562d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
15639566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitSetRuntimeSplits_Private(pc));
15649566063dSJacob Faibussowitsch     if (jac->splitdefined) PetscCall(PetscInfo(pc, "Splits defined using the options database\n"));
1565d32f9abdSBarry Smith   }
1566c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
15679566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetEnum(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_fieldsplit_schur_factorization_type", PCFieldSplitSchurFactTypes, (PetscEnum *)&jac->schurfactorization, &flg));
15689566063dSJacob Faibussowitsch     if (flg) PetscCall(PetscInfo(pc, "Deprecated use of -pc_fieldsplit_schur_factorization_type\n"));
15699566063dSJacob Faibussowitsch     PetscCall(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));
15709566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-pc_fieldsplit_schur_precondition", "How to build preconditioner for Schur complement", "PCFieldSplitSetSchurPre", PCFieldSplitSchurPreTypes, (PetscEnum)jac->schurpre, (PetscEnum *)&jac->schurpre, NULL));
15719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-pc_fieldsplit_schur_scale", "Scale Schur complement", "PCFieldSplitSetSchurScale", jac->schurscale, &jac->schurscale, NULL));
1572a51937d4SCarola Kruse   } else if (jac->type == PC_COMPOSITE_GKB) {
15739566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_fieldsplit_gkb_tol", "The tolerance for the lower bound stopping criterion", "PCFieldSplitGKBTol", jac->gkbtol, &jac->gkbtol, NULL));
15749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_fieldsplit_gkb_delay", "The delay value for lower bound criterion", "PCFieldSplitGKBDelay", jac->gkbdelay, &jac->gkbdelay, NULL));
15759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_fieldsplit_gkb_nu", "Parameter in augmented Lagrangian approach", "PCFieldSplitGKBNu", jac->gkbnu, &jac->gkbnu, NULL));
157663a3b9bcSJacob Faibussowitsch     PetscCheck(jac->gkbnu >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nu cannot be less than 0: value %g", (double)jac->gkbnu);
15779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_fieldsplit_gkb_maxit", "Maximum allowed number of iterations", "PCFieldSplitGKBMaxit", jac->gkbmaxit, &jac->gkbmaxit, NULL));
15789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_fieldsplit_gkb_monitor", "Prints number of GKB iterations and error", "PCFieldSplitGKB", jac->gkbmonitor, &jac->gkbmonitor, NULL));
1579c5d2311dSJed Brown   }
158034603f55SBarry Smith   /*
158134603f55SBarry Smith     In the initial call to this routine the sub-solver data structures do not exist so we cannot call KSPSetFromOptions() on them yet.
158234603f55SBarry Smith     But after the initial setup of ALL the layers of sub-solvers is completed we do want to call KSPSetFromOptions() on the sub-solvers every time it
158334603f55SBarry Smith     is called on the outer solver in case changes were made in the options database
158434603f55SBarry Smith 
158534603f55SBarry Smith     But even after PCSetUp_FieldSplit() is called all the options inside the inner levels of sub-solvers may still not have been set thus we only call the KSPSetFromOptions()
158634603f55SBarry Smith     if we know that the entire stack of sub-solvers below this have been complete instantiated, we check this by seeing if any solver iterations are complete.
158734603f55SBarry Smith     Without this extra check test p2p1fetidp_olof_full and others fail with incorrect matrix types.
158834603f55SBarry Smith 
158934603f55SBarry Smith     There could be a negative side effect of calling the KSPSetFromOptions() below.
159034603f55SBarry Smith 
159134603f55SBarry Smith     If one captured the PetscObjectState of the options database one could skip these calls if the database has not changed from the previous call
159234603f55SBarry Smith   */
159334603f55SBarry Smith   if (jac->issetup) {
159434603f55SBarry Smith     PC_FieldSplitLink ilink = jac->head;
159534603f55SBarry Smith     if (jac->type == PC_COMPOSITE_SCHUR) {
159634603f55SBarry Smith       if (jac->kspupper && jac->kspupper->totalits > 0) PetscCall(KSPSetFromOptions(jac->kspupper));
159734603f55SBarry Smith       if (jac->kspschur && jac->kspschur->totalits > 0) PetscCall(KSPSetFromOptions(jac->kspschur));
159834603f55SBarry Smith     }
159934603f55SBarry Smith     while (ilink) {
160034603f55SBarry Smith       if (ilink->ksp->totalits > 0) PetscCall(KSPSetFromOptions(ilink->ksp));
160134603f55SBarry Smith       ilink = ilink->next;
160234603f55SBarry Smith     }
160334603f55SBarry Smith   }
1604d0609cedSBarry Smith   PetscOptionsHeadEnd();
16050971522cSBarry Smith   PetscFunctionReturn(0);
16060971522cSBarry Smith }
16070971522cSBarry Smith 
16089371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc, const char splitname[], PetscInt n, const PetscInt *fields, const PetscInt *fields_col) {
160997bbdb24SBarry Smith   PC_FieldSplit    *jac = (PC_FieldSplit *)pc->data;
16105a9f2f41SSatish Balay   PC_FieldSplitLink ilink, next = jac->head;
161169a612a9SBarry Smith   char              prefix[128];
16125d4c12cdSJungho Lee   PetscInt          i;
16130971522cSBarry Smith 
16140971522cSBarry Smith   PetscFunctionBegin;
16156c924f48SJed Brown   if (jac->splitdefined) {
16169566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Ignoring new split \"%s\" because the splits have already been defined\n", splitname));
16176c924f48SJed Brown     PetscFunctionReturn(0);
16186c924f48SJed Brown   }
161951f519a2SBarry Smith   for (i = 0; i < n; i++) {
162063a3b9bcSJacob Faibussowitsch     PetscCheck(fields[i] < jac->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " requested but only %" PetscInt_FMT " exist", fields[i], jac->bs);
162163a3b9bcSJacob Faibussowitsch     PetscCheck(fields[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %" PetscInt_FMT " requested", fields[i]);
162251f519a2SBarry Smith   }
16239566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
1624a04f6461SBarry Smith   if (splitname) {
16259566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(splitname, &ilink->splitname));
1626a04f6461SBarry Smith   } else {
16279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(3, &ilink->splitname));
162863a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ilink->splitname, 2, "%" PetscInt_FMT, jac->nsplits));
1629a04f6461SBarry Smith   }
16309df09d43SBarry 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 */
16319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &ilink->fields));
16329566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ilink->fields, fields, n));
16339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &ilink->fields_col));
16349566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ilink->fields_col, fields_col, n));
16352fa5cd67SKarl Rupp 
16365a9f2f41SSatish Balay   ilink->nfields = n;
16370298fd71SBarry Smith   ilink->next    = NULL;
16389566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &ilink->ksp));
16399566063dSJacob Faibussowitsch   PetscCall(KSPSetErrorIfNotConverged(ilink->ksp, pc->erroriffailure));
16409566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->ksp, (PetscObject)pc, 1));
16419566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ilink->ksp, KSPPREONLY));
164269a612a9SBarry Smith 
16439566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname));
16449566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(ilink->ksp, prefix));
16450971522cSBarry Smith 
16460971522cSBarry Smith   if (!next) {
16475a9f2f41SSatish Balay     jac->head       = ilink;
16480298fd71SBarry Smith     ilink->previous = NULL;
16490971522cSBarry Smith   } else {
1650ad540459SPierre Jolivet     while (next->next) next = next->next;
16515a9f2f41SSatish Balay     next->next      = ilink;
165251f519a2SBarry Smith     ilink->previous = next;
16530971522cSBarry Smith   }
16540971522cSBarry Smith   jac->nsplits++;
16550971522cSBarry Smith   PetscFunctionReturn(0);
16560971522cSBarry Smith }
16570971522cSBarry Smith 
16589371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSchurGetSubKSP_FieldSplit(PC pc, PetscInt *n, KSP **subksp) {
1659285fb4e2SStefano Zampini   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
1660285fb4e2SStefano Zampini 
1661285fb4e2SStefano Zampini   PetscFunctionBegin;
1662285fb4e2SStefano Zampini   *subksp = NULL;
1663285fb4e2SStefano Zampini   if (n) *n = 0;
1664285fb4e2SStefano Zampini   if (jac->type == PC_COMPOSITE_SCHUR) {
1665285fb4e2SStefano Zampini     PetscInt nn;
1666285fb4e2SStefano Zampini 
166728b400f6SJacob Faibussowitsch     PetscCheck(jac->schur, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitSchurGetSubKSP()");
166863a3b9bcSJacob Faibussowitsch     PetscCheck(jac->nsplits == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Unexpected number of splits %" PetscInt_FMT " != 2", jac->nsplits);
1669285fb4e2SStefano Zampini     nn = jac->nsplits + (jac->kspupper != jac->head->ksp ? 1 : 0);
16709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nn, subksp));
1671285fb4e2SStefano Zampini     (*subksp)[0] = jac->head->ksp;
1672285fb4e2SStefano Zampini     (*subksp)[1] = jac->kspschur;
1673285fb4e2SStefano Zampini     if (jac->kspupper != jac->head->ksp) (*subksp)[2] = jac->kspupper;
1674285fb4e2SStefano Zampini     if (n) *n = nn;
1675285fb4e2SStefano Zampini   }
1676285fb4e2SStefano Zampini   PetscFunctionReturn(0);
1677285fb4e2SStefano Zampini }
1678285fb4e2SStefano Zampini 
16799371c9d4SSatish Balay static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc, PetscInt *n, KSP **subksp) {
1680e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
1681e69d4d44SBarry Smith 
1682e69d4d44SBarry Smith   PetscFunctionBegin;
168328b400f6SJacob Faibussowitsch   PetscCheck(jac->schur, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()");
16849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jac->nsplits, subksp));
16859566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetKSP(jac->schur, *subksp));
16862fa5cd67SKarl Rupp 
1687e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
168813e0d083SBarry Smith   if (n) *n = jac->nsplits;
1689e69d4d44SBarry Smith   PetscFunctionReturn(0);
1690e69d4d44SBarry Smith }
16910971522cSBarry Smith 
16929371c9d4SSatish Balay static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc, PetscInt *n, KSP **subksp) {
16930971522cSBarry Smith   PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
16940971522cSBarry Smith   PetscInt          cnt   = 0;
16955a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
16960971522cSBarry Smith 
16970971522cSBarry Smith   PetscFunctionBegin;
16989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jac->nsplits, subksp));
16995a9f2f41SSatish Balay   while (ilink) {
17005a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
17015a9f2f41SSatish Balay     ilink            = ilink->next;
17020971522cSBarry Smith   }
170363a3b9bcSJacob Faibussowitsch   PetscCheck(cnt == jac->nsplits, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Corrupt PCFIELDSPLIT object: number of splits in linked list %" PetscInt_FMT " does not match number in object %" PetscInt_FMT, cnt, jac->nsplits);
170413e0d083SBarry Smith   if (n) *n = jac->nsplits;
17050971522cSBarry Smith   PetscFunctionReturn(0);
17060971522cSBarry Smith }
17070971522cSBarry Smith 
17086dbb499eSCian Wilson /*@C
1709f1580f4eSBarry Smith     PCFieldSplitRestrictIS - Restricts the fieldsplit `IS`s to be within a given `IS`.
17106dbb499eSCian Wilson 
17116dbb499eSCian Wilson     Input Parameters:
17126dbb499eSCian Wilson +   pc  - the preconditioner context
1713a2b725a8SWilliam Gropp -   is - the index set that defines the indices to which the fieldsplit is to be restricted
17146dbb499eSCian Wilson 
17156dbb499eSCian Wilson     Level: advanced
17166dbb499eSCian Wilson 
1717f1580f4eSBarry Smith     Developer Note:
1718f1580f4eSBarry Smith     It seems the resulting `IS`s will not cover the entire space, so
1719f1580f4eSBarry Smith     how can they define a convergent preconditioner? Needs explaining.
1720f1580f4eSBarry Smith 
1721f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()`
17226dbb499eSCian Wilson @*/
17239371c9d4SSatish Balay PetscErrorCode PCFieldSplitRestrictIS(PC pc, IS isy) {
17246dbb499eSCian Wilson   PetscFunctionBegin;
17256dbb499eSCian Wilson   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
17266dbb499eSCian Wilson   PetscValidHeaderSpecific(isy, IS_CLASSID, 2);
1727cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitRestrictIS_C", (PC, IS), (pc, isy));
17286dbb499eSCian Wilson   PetscFunctionReturn(0);
17296dbb499eSCian Wilson }
17306dbb499eSCian Wilson 
17319371c9d4SSatish Balay static PetscErrorCode PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy) {
17326dbb499eSCian Wilson   PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
17336dbb499eSCian Wilson   PC_FieldSplitLink ilink = jac->head, next;
17346dbb499eSCian Wilson   PetscInt          localsize, size, sizez, i;
17356dbb499eSCian Wilson   const PetscInt   *ind, *indz;
17366dbb499eSCian Wilson   PetscInt         *indc, *indcz;
17376dbb499eSCian Wilson   PetscBool         flg;
17386dbb499eSCian Wilson 
17396dbb499eSCian Wilson   PetscFunctionBegin;
17409566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isy, &localsize));
17419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Scan(&localsize, &size, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)isy)));
17426dbb499eSCian Wilson   size -= localsize;
17436dbb499eSCian Wilson   while (ilink) {
17446dbb499eSCian Wilson     IS isrl, isr;
17451c7cfba8SBarry Smith     PC subpc;
17469566063dSJacob Faibussowitsch     PetscCall(ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl));
17479566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isrl, &localsize));
17489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(localsize, &indc));
17499566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(isrl, &ind));
17509566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(indc, ind, localsize));
17519566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(isrl, &ind));
17529566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isrl));
17536dbb499eSCian Wilson     for (i = 0; i < localsize; i++) *(indc + i) += size;
17549566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)isy), localsize, indc, PETSC_OWN_POINTER, &isr));
17559566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)isr));
17569566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ilink->is));
17576dbb499eSCian Wilson     ilink->is = isr;
17589566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)isr));
17599566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ilink->is_col));
17606dbb499eSCian Wilson     ilink->is_col = isr;
17619566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isr));
17629566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ilink->ksp, &subpc));
17639566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)subpc, PCFIELDSPLIT, &flg));
17646dbb499eSCian Wilson     if (flg) {
17656dbb499eSCian Wilson       IS       iszl, isz;
17666dbb499eSCian Wilson       MPI_Comm comm;
17679566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(ilink->is, &localsize));
17686dbb499eSCian Wilson       comm = PetscObjectComm((PetscObject)ilink->is);
17699566063dSJacob Faibussowitsch       PetscCall(ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl));
17709566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&localsize, &sizez, 1, MPIU_INT, MPI_SUM, comm));
17716dbb499eSCian Wilson       sizez -= localsize;
17729566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(iszl, &localsize));
17739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(localsize, &indcz));
17749566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(iszl, &indz));
17759566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(indcz, indz, localsize));
17769566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(iszl, &indz));
17779566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&iszl));
17786dbb499eSCian Wilson       for (i = 0; i < localsize; i++) *(indcz + i) += sizez;
17799566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, localsize, indcz, PETSC_OWN_POINTER, &isz));
17809566063dSJacob Faibussowitsch       PetscCall(PCFieldSplitRestrictIS(subpc, isz));
17819566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isz));
17826dbb499eSCian Wilson     }
17836dbb499eSCian Wilson     next  = ilink->next;
17846dbb499eSCian Wilson     ilink = next;
17856dbb499eSCian Wilson   }
17866dbb499eSCian Wilson   jac->isrestrict = PETSC_TRUE;
17876dbb499eSCian Wilson   PetscFunctionReturn(0);
17886dbb499eSCian Wilson }
17896dbb499eSCian Wilson 
17909371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc, const char splitname[], IS is) {
1791704ba839SBarry Smith   PC_FieldSplit    *jac = (PC_FieldSplit *)pc->data;
1792704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1793704ba839SBarry Smith   char              prefix[128];
1794704ba839SBarry Smith 
1795704ba839SBarry Smith   PetscFunctionBegin;
17966c924f48SJed Brown   if (jac->splitdefined) {
17979566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Ignoring new split \"%s\" because the splits have already been defined\n", splitname));
17986c924f48SJed Brown     PetscFunctionReturn(0);
17996c924f48SJed Brown   }
18009566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
1801a04f6461SBarry Smith   if (splitname) {
18029566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(splitname, &ilink->splitname));
1803a04f6461SBarry Smith   } else {
18049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(8, &ilink->splitname));
180563a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ilink->splitname, 7, "%" PetscInt_FMT, jac->nsplits));
1806a04f6461SBarry Smith   }
18079df09d43SBarry 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 */
18089566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)is));
18099566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ilink->is));
1810b5787286SJed Brown   ilink->is = is;
18119566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)is));
18129566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ilink->is_col));
1813b5787286SJed Brown   ilink->is_col = is;
18140298fd71SBarry Smith   ilink->next   = NULL;
18159566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &ilink->ksp));
18169566063dSJacob Faibussowitsch   PetscCall(KSPSetErrorIfNotConverged(ilink->ksp, pc->erroriffailure));
18179566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->ksp, (PetscObject)pc, 1));
18189566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ilink->ksp, KSPPREONLY));
1819704ba839SBarry Smith 
18209566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname));
18219566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(ilink->ksp, prefix));
1822704ba839SBarry Smith 
1823704ba839SBarry Smith   if (!next) {
1824704ba839SBarry Smith     jac->head       = ilink;
18250298fd71SBarry Smith     ilink->previous = NULL;
1826704ba839SBarry Smith   } else {
1827ad540459SPierre Jolivet     while (next->next) next = next->next;
1828704ba839SBarry Smith     next->next      = ilink;
1829704ba839SBarry Smith     ilink->previous = next;
1830704ba839SBarry Smith   }
1831704ba839SBarry Smith   jac->nsplits++;
1832704ba839SBarry Smith   PetscFunctionReturn(0);
1833704ba839SBarry Smith }
1834704ba839SBarry Smith 
183594ef8ddeSSatish Balay /*@C
1836f1580f4eSBarry Smith     PCFieldSplitSetFields - Sets the fields that define one particular split in the field split preconditioner
18370971522cSBarry Smith 
1838f1580f4eSBarry Smith     Logically Collective on pc
18390971522cSBarry Smith 
18400971522cSBarry Smith     Input Parameters:
18410971522cSBarry Smith +   pc  - the preconditioner context
18420298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
18430971522cSBarry Smith .   n - the number of fields in this split
1844db4c96c1SJed Brown -   fields - the fields in this split
18450971522cSBarry Smith 
18460971522cSBarry Smith     Level: intermediate
18470971522cSBarry Smith 
184895452b02SPatrick Sanan     Notes:
1849f1580f4eSBarry Smith     Use `PCFieldSplitSetIS()` to set a  general set of indices as a split.
1850d32f9abdSBarry Smith 
1851f1580f4eSBarry Smith      `PCFieldSplitSetFields()` is for defining fields as strided blocks. For example, if the block
1852f1580f4eSBarry Smith      size is three then one can define a split as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean
1853d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1854f1580f4eSBarry Smith      where the numbered entries indicate what is in the split.
1855d32f9abdSBarry Smith 
1856db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1857db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1858db4c96c1SJed Brown 
1859f1580f4eSBarry Smith    Developer Note:
1860f1580f4eSBarry Smith    This routine does not actually create the `IS` representing the split, that is delayed
1861f1580f4eSBarry Smith    until `PCSetUp_FieldSplit()`, because information about the vector/matrix layouts may not be
18625d4c12cdSJungho Lee    available when this routine is called.
18635d4c12cdSJungho Lee 
1864f1580f4eSBarry Smith .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetBlockSize()`, `PCFieldSplitSetIS()`, `PCFieldSplitRestrictIS()`
18650971522cSBarry Smith @*/
18669371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetFields(PC pc, const char splitname[], PetscInt n, const PetscInt *fields, const PetscInt *fields_col) {
18670971522cSBarry Smith   PetscFunctionBegin;
18680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1869db4c96c1SJed Brown   PetscValidCharPointer(splitname, 2);
187063a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Provided number of fields %" PetscInt_FMT " in split \"%s\" not positive", n, splitname);
1871064a246eSJacob Faibussowitsch   PetscValidIntPointer(fields, 4);
1872cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetFields_C", (PC, const char[], PetscInt, const PetscInt *, const PetscInt *), (pc, splitname, n, fields, fields_col));
18730971522cSBarry Smith   PetscFunctionReturn(0);
18740971522cSBarry Smith }
18750971522cSBarry Smith 
1876c84da90fSDmitry Karpeev /*@
1877f1580f4eSBarry Smith     PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat) to build
1878f1580f4eSBarry Smith     the sub-matrices associated with each split. Where `KSPSetOperators`(ksp,Amat,Pmat)) was used to supply the operators.
1879c84da90fSDmitry Karpeev 
1880f1580f4eSBarry Smith     Logically Collective on pc
1881c84da90fSDmitry Karpeev 
1882c84da90fSDmitry Karpeev     Input Parameters:
1883c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1884c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1885c84da90fSDmitry Karpeev 
1886f1580f4eSBarry Smith     Options Database Keys:
1887147403d9SBarry Smith .   -pc_fieldsplit_diag_use_amat - use the Amat to provide the diagonal blocks
1888c84da90fSDmitry Karpeev 
1889c84da90fSDmitry Karpeev     Level: intermediate
1890c84da90fSDmitry Karpeev 
1891f1580f4eSBarry Smith .seealso: `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitGetDiagUseAmat()`, `PCFieldSplitSetOffDiagUseAmat()`, `PCFIELDSPLIT`
1892c84da90fSDmitry Karpeev @*/
18939371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetDiagUseAmat(PC pc, PetscBool flg) {
1894c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
1895c84da90fSDmitry Karpeev   PetscBool      isfs;
1896c84da90fSDmitry Karpeev 
1897c84da90fSDmitry Karpeev   PetscFunctionBegin;
1898c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
18999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs));
190028b400f6SJacob Faibussowitsch   PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT);
1901c84da90fSDmitry Karpeev   jac->diag_use_amat = flg;
1902c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1903c84da90fSDmitry Karpeev }
1904c84da90fSDmitry Karpeev 
1905c84da90fSDmitry Karpeev /*@
1906f1580f4eSBarry Smith     PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat) to build
1907f1580f4eSBarry Smith     the sub-matrices associated with each split.  Where `KSPSetOperators`(ksp,Amat,Pmat)) was used to supply the operators.
1908c84da90fSDmitry Karpeev 
1909f1580f4eSBarry Smith     Logically Collective on pc
1910c84da90fSDmitry Karpeev 
1911c84da90fSDmitry Karpeev     Input Parameters:
1912c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1913c84da90fSDmitry Karpeev 
1914c84da90fSDmitry Karpeev     Output Parameters:
1915c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1916c84da90fSDmitry Karpeev 
1917c84da90fSDmitry Karpeev     Level: intermediate
1918c84da90fSDmitry Karpeev 
1919f1580f4eSBarry Smith .seealso: `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitSetDiagUseAmat()`, `PCFieldSplitGetOffDiagUseAmat()`, `PCFIELDSPLIT`
1920c84da90fSDmitry Karpeev @*/
19219371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetDiagUseAmat(PC pc, PetscBool *flg) {
1922c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
1923c84da90fSDmitry Karpeev   PetscBool      isfs;
1924c84da90fSDmitry Karpeev 
1925c84da90fSDmitry Karpeev   PetscFunctionBegin;
1926c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1927534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
19289566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs));
192928b400f6SJacob Faibussowitsch   PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT);
1930c84da90fSDmitry Karpeev   *flg = jac->diag_use_amat;
1931c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1932c84da90fSDmitry Karpeev }
1933c84da90fSDmitry Karpeev 
1934c84da90fSDmitry Karpeev /*@
1935f1580f4eSBarry Smith     PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat) to build
1936f1580f4eSBarry Smith     the sub-matrices associated with each split.  Where `KSPSetOperators`(ksp,Amat,Pmat)) was used to supply the operators.
1937c84da90fSDmitry Karpeev 
1938f1580f4eSBarry Smith     Logically Collective on pc
1939c84da90fSDmitry Karpeev 
1940c84da90fSDmitry Karpeev     Input Parameters:
1941c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1942c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1943c84da90fSDmitry Karpeev 
1944f1580f4eSBarry Smith     Options Database Keys:
1945147403d9SBarry Smith .     -pc_fieldsplit_off_diag_use_amat <bool> - use the Amat to extract the off-diagonal blocks
1946c84da90fSDmitry Karpeev 
1947c84da90fSDmitry Karpeev     Level: intermediate
1948c84da90fSDmitry Karpeev 
1949f1580f4eSBarry Smith .seealso: `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitGetOffDiagUseAmat()`, `PCFieldSplitSetDiagUseAmat()`, `PCFIELDSPLIT`
1950c84da90fSDmitry Karpeev @*/
19519371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC pc, PetscBool flg) {
1952c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
1953c84da90fSDmitry Karpeev   PetscBool      isfs;
1954c84da90fSDmitry Karpeev 
1955c84da90fSDmitry Karpeev   PetscFunctionBegin;
1956c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs));
195828b400f6SJacob Faibussowitsch   PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT);
1959c84da90fSDmitry Karpeev   jac->offdiag_use_amat = flg;
1960c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1961c84da90fSDmitry Karpeev }
1962c84da90fSDmitry Karpeev 
1963c84da90fSDmitry Karpeev /*@
1964f1580f4eSBarry Smith     PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat) to build
1965f1580f4eSBarry Smith     the sub-matrices associated with each split.  Where `KSPSetOperators`(ksp,Amat,Pmat)) was used to supply the operators.
1966c84da90fSDmitry Karpeev 
1967f1580f4eSBarry Smith     Logically Collective on pc
1968c84da90fSDmitry Karpeev 
1969c84da90fSDmitry Karpeev     Input Parameters:
1970c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1971c84da90fSDmitry Karpeev 
1972c84da90fSDmitry Karpeev     Output Parameters:
1973c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1974c84da90fSDmitry Karpeev 
1975c84da90fSDmitry Karpeev     Level: intermediate
1976c84da90fSDmitry Karpeev 
1977f1580f4eSBarry Smith .seealso: `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitSetOffDiagUseAmat()`, `PCFieldSplitGetDiagUseAmat()`, `PCFIELDSPLIT`
1978c84da90fSDmitry Karpeev @*/
19799371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC pc, PetscBool *flg) {
1980c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
1981c84da90fSDmitry Karpeev   PetscBool      isfs;
1982c84da90fSDmitry Karpeev 
1983c84da90fSDmitry Karpeev   PetscFunctionBegin;
1984c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1985534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
19869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs));
198728b400f6SJacob Faibussowitsch   PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT);
1988c84da90fSDmitry Karpeev   *flg = jac->offdiag_use_amat;
1989c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1990c84da90fSDmitry Karpeev }
1991c84da90fSDmitry Karpeev 
1992218d4943SBarry Smith /*@C
1993f1580f4eSBarry Smith     PCFieldSplitSetIS - Sets the exact elements for a split in a `PCFIELDSPLIT`
1994704ba839SBarry Smith 
1995f1580f4eSBarry Smith     Logically Collective on pc
1996704ba839SBarry Smith 
1997704ba839SBarry Smith     Input Parameters:
1998704ba839SBarry Smith +   pc  - the preconditioner context
19990298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
2000f1580f4eSBarry Smith -   is - the index set that defines the elements in this split
2001704ba839SBarry Smith 
2002a6ffb8dbSJed Brown     Notes:
2003f1580f4eSBarry Smith     Use `PCFieldSplitSetFields()`, for splits defined by strided types.
2004a6ffb8dbSJed Brown 
2005db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
2006db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
2007d32f9abdSBarry Smith 
2008704ba839SBarry Smith     Level: intermediate
2009704ba839SBarry Smith 
2010db781477SPatrick Sanan .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetBlockSize()`
2011704ba839SBarry Smith @*/
20129371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetIS(PC pc, const char splitname[], IS is) {
2013704ba839SBarry Smith   PetscFunctionBegin;
20140700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20157b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname, 2);
2016db4c96c1SJed Brown   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
2017cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetIS_C", (PC, const char[], IS), (pc, splitname, is));
2018704ba839SBarry Smith   PetscFunctionReturn(0);
2019704ba839SBarry Smith }
2020704ba839SBarry Smith 
202194ef8ddeSSatish Balay /*@C
2022f1580f4eSBarry Smith     PCFieldSplitGetIS - Retrieves the elements for a split as an `IS`
202357a9adfeSMatthew G Knepley 
2024f1580f4eSBarry Smith     Logically Collective on pc
202557a9adfeSMatthew G Knepley 
202657a9adfeSMatthew G Knepley     Input Parameters:
202757a9adfeSMatthew G Knepley +   pc  - the preconditioner context
202857a9adfeSMatthew G Knepley -   splitname - name of this split
202957a9adfeSMatthew G Knepley 
203057a9adfeSMatthew G Knepley     Output Parameter:
2031f1580f4eSBarry Smith -   is - the index set that defines the elements in this split, or NULL if the split is not found
203257a9adfeSMatthew G Knepley 
203357a9adfeSMatthew G Knepley     Level: intermediate
203457a9adfeSMatthew G Knepley 
2035db781477SPatrick Sanan .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetIS()`
203657a9adfeSMatthew G Knepley @*/
20379371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetIS(PC pc, const char splitname[], IS *is) {
203857a9adfeSMatthew G Knepley   PetscFunctionBegin;
203957a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
204057a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname, 2);
204157a9adfeSMatthew G Knepley   PetscValidPointer(is, 3);
204257a9adfeSMatthew G Knepley   {
204357a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
204457a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
204557a9adfeSMatthew G Knepley     PetscBool         found;
204657a9adfeSMatthew G Knepley 
20470298fd71SBarry Smith     *is = NULL;
204857a9adfeSMatthew G Knepley     while (ilink) {
20499566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(ilink->splitname, splitname, &found));
205057a9adfeSMatthew G Knepley       if (found) {
205157a9adfeSMatthew G Knepley         *is = ilink->is;
205257a9adfeSMatthew G Knepley         break;
205357a9adfeSMatthew G Knepley       }
205457a9adfeSMatthew G Knepley       ilink = ilink->next;
205557a9adfeSMatthew G Knepley     }
205657a9adfeSMatthew G Knepley   }
205757a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
205857a9adfeSMatthew G Knepley }
205957a9adfeSMatthew G Knepley 
2060998f007dSPierre Jolivet /*@C
2061f1580f4eSBarry Smith     PCFieldSplitGetISByIndex - Retrieves the elements for a given split as an `IS`
2062998f007dSPierre Jolivet 
2063f1580f4eSBarry Smith     Logically Collective on pc
2064998f007dSPierre Jolivet 
2065998f007dSPierre Jolivet     Input Parameters:
2066998f007dSPierre Jolivet +   pc  - the preconditioner context
2067998f007dSPierre Jolivet -   index - index of this split
2068998f007dSPierre Jolivet 
2069998f007dSPierre Jolivet     Output Parameter:
2070f1580f4eSBarry Smith -   is - the index set that defines the elements in this split
2071998f007dSPierre Jolivet 
2072998f007dSPierre Jolivet     Level: intermediate
2073998f007dSPierre Jolivet 
2074db781477SPatrick Sanan .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitGetIS()`, `PCFieldSplitSetIS()`
2075998f007dSPierre Jolivet @*/
20769371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetISByIndex(PC pc, PetscInt index, IS *is) {
2077998f007dSPierre Jolivet   PetscFunctionBegin;
207863a3b9bcSJacob Faibussowitsch   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %" PetscInt_FMT " requested", index);
2079998f007dSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2080998f007dSPierre Jolivet   PetscValidPointer(is, 3);
2081998f007dSPierre Jolivet   {
2082998f007dSPierre Jolivet     PC_FieldSplit    *jac   = (PC_FieldSplit *)pc->data;
2083998f007dSPierre Jolivet     PC_FieldSplitLink ilink = jac->head;
2084998f007dSPierre Jolivet     PetscInt          i     = 0;
208563a3b9bcSJacob Faibussowitsch     PetscCheck(index < jac->nsplits, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " requested but only %" PetscInt_FMT " exist", index, jac->nsplits);
2086998f007dSPierre Jolivet 
2087998f007dSPierre Jolivet     while (i < index) {
2088998f007dSPierre Jolivet       ilink = ilink->next;
2089998f007dSPierre Jolivet       ++i;
2090998f007dSPierre Jolivet     }
20919566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitGetIS(pc, ilink->splitname, is));
2092998f007dSPierre Jolivet   }
2093998f007dSPierre Jolivet   PetscFunctionReturn(0);
2094998f007dSPierre Jolivet }
2095998f007dSPierre Jolivet 
209651f519a2SBarry Smith /*@
209751f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
2098f1580f4eSBarry Smith       fieldsplit preconditioner when calling `PCFieldSplitSetIS()`. If not set the matrix block size is used.
209951f519a2SBarry Smith 
2100f1580f4eSBarry Smith     Logically Collective on pc
210151f519a2SBarry Smith 
210251f519a2SBarry Smith     Input Parameters:
210351f519a2SBarry Smith +   pc  - the preconditioner context
210451f519a2SBarry Smith -   bs - the block size
210551f519a2SBarry Smith 
210651f519a2SBarry Smith     Level: intermediate
210751f519a2SBarry Smith 
2108f1580f4eSBarry Smith .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()`
210951f519a2SBarry Smith @*/
21109371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetBlockSize(PC pc, PetscInt bs) {
211151f519a2SBarry Smith   PetscFunctionBegin;
21120700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2113c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, bs, 2);
2114cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetBlockSize_C", (PC, PetscInt), (pc, bs));
211551f519a2SBarry Smith   PetscFunctionReturn(0);
211651f519a2SBarry Smith }
211751f519a2SBarry Smith 
21180971522cSBarry Smith /*@C
2119f1580f4eSBarry Smith    PCFieldSplitGetSubKSP - Gets the `KSP` contexts for all splits
21200971522cSBarry Smith 
2121f1580f4eSBarry Smith    Collective on pc
21220971522cSBarry Smith 
21230971522cSBarry Smith    Input Parameter:
21240971522cSBarry Smith .  pc - the preconditioner context
21250971522cSBarry Smith 
21260971522cSBarry Smith    Output Parameters:
212713e0d083SBarry Smith +  n - the number of splits
2128f1580f4eSBarry Smith -  subksp - the array of `KSP` contexts
2129196cc216SBarry Smith 
21300971522cSBarry Smith    Level: advanced
21310971522cSBarry Smith 
2132f1580f4eSBarry Smith    Notes:
2133f1580f4eSBarry Smith    After `PCFieldSplitGetSubKSP()` the array of `KSP`s is to be freed by the user with `PetscFree()`
2134f1580f4eSBarry Smith    (not the `KSP`, just the array that contains them).
2135f1580f4eSBarry Smith 
2136f1580f4eSBarry Smith    You must call `PCSetUp()` before calling `PCFieldSplitGetSubKSP()`.
2137f1580f4eSBarry Smith 
2138f1580f4eSBarry Smith    If the fieldsplit is of type `PC_COMPOSITE_SCHUR`, it returns the `KSP` object used inside the
2139f1580f4eSBarry Smith    Schur complement and the `KSP` object used to iterate over the Schur complement.
2140f1580f4eSBarry Smith    To access all the `KSP` objects used in `PC_COMPOSITE_SCHUR`, use `PCFieldSplitSchurGetSubKSP()`.
2141f1580f4eSBarry Smith 
2142f1580f4eSBarry Smith    If the fieldsplit is of type `PC_COMPOSITE_GKB`, it returns the `KSP` object used to solve the
2143f1580f4eSBarry Smith    inner linear system defined by the matrix H in each loop.
2144f1580f4eSBarry Smith 
2145f1580f4eSBarry Smith    Fortran Usage:
2146f1580f4eSBarry Smith    You must pass in a `KSP` array that is large enough to contain all the `KSP`s.
2147f1580f4eSBarry Smith    You can call `PCFieldSplitGetSubKSP`(pc,n,`PETSC_NULL_KSP`,ierr) to determine how large the
2148f1580f4eSBarry Smith    `KSP` array must be.
2149f1580f4eSBarry Smith 
2150f1580f4eSBarry Smith    Developer Note:
2151f1580f4eSBarry Smith    There should be a `PCFieldSplitRestoreSubKSP()` instead of requiring the user to call `PetscFree()`
2152f1580f4eSBarry Smith 
2153f1580f4eSBarry Smith    The Fortran interface should be modernized to return directly the array of values.
2154f1580f4eSBarry Smith 
2155f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()`, `PCFieldSplitSchurGetSubKSP()`
21560971522cSBarry Smith @*/
21579371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetSubKSP(PC pc, PetscInt *n, KSP *subksp[]) {
21580971522cSBarry Smith   PetscFunctionBegin;
21590700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
216013e0d083SBarry Smith   if (n) PetscValidIntPointer(n, 2);
2161cac4c232SBarry Smith   PetscUseMethod(pc, "PCFieldSplitGetSubKSP_C", (PC, PetscInt *, KSP **), (pc, n, subksp));
21620971522cSBarry Smith   PetscFunctionReturn(0);
21630971522cSBarry Smith }
21640971522cSBarry Smith 
2165285fb4e2SStefano Zampini /*@C
2166f1580f4eSBarry Smith    PCFieldSplitSchurGetSubKSP - Gets the `KSP` contexts used inside the Schur complement based `PCFIELDSPLIT`
2167285fb4e2SStefano Zampini 
2168f1580f4eSBarry Smith    Collective on `KSP`
2169285fb4e2SStefano Zampini 
2170285fb4e2SStefano Zampini    Input Parameter:
2171285fb4e2SStefano Zampini .  pc - the preconditioner context
2172285fb4e2SStefano Zampini 
2173285fb4e2SStefano Zampini    Output Parameters:
2174285fb4e2SStefano Zampini +  n - the number of splits
2175f1580f4eSBarry Smith -  subksp - the array of `KSP` contexts
2176285fb4e2SStefano Zampini 
2177285fb4e2SStefano Zampini    Level: advanced
2178285fb4e2SStefano Zampini 
2179f1580f4eSBarry Smith    Notes:
2180f1580f4eSBarry Smith    After `PCFieldSplitSchurGetSubKSP()` the array of `KSP`s is to be freed by the user with `PetscFree()`
2181f1580f4eSBarry Smith    (not the `KSP` just the array that contains them).
2182f1580f4eSBarry Smith 
2183f1580f4eSBarry Smith    You must call `PCSetUp()` before calling `PCFieldSplitSchurGetSubKSP()`.
2184f1580f4eSBarry Smith 
2185f1580f4eSBarry Smith    If the fieldsplit type is of type `PC_COMPOSITE_SCHUR`, it returns (in order)
2186f1580f4eSBarry Smith +  1  - the `KSP` used for the (1,1) block
2187f1580f4eSBarry Smith .  2  - the `KSP` used for the Schur complement (not the one used for the interior Schur solver)
2188f1580f4eSBarry Smith -  3  - the `KSP` used for the (1,1) block in the upper triangular factor (if different from that of the (1,1) block).
2189f1580f4eSBarry Smith 
2190f1580f4eSBarry Smith    It returns a null array if the fieldsplit is not of type `PC_COMPOSITE_SCHUR`; in this case, you should use `PCFieldSplitGetSubKSP()`.
2191f1580f4eSBarry Smith 
2192f1580f4eSBarry Smith    Fortran Note:
2193f1580f4eSBarry Smith    You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
2194f1580f4eSBarry Smith    You can call `PCFieldSplitSchurGetSubKSP`(pc,n,`PETSC_NULL_KSP`,ierr) to determine how large the
2195f1580f4eSBarry Smith    `KSP` array must be.
2196f1580f4eSBarry Smith 
2197f1580f4eSBarry Smith    Developer Notes:
2198f1580f4eSBarry Smith    There should be a `PCFieldSplitRestoreSubKSP()` instead of requiring the user to call `PetscFree()`
2199f1580f4eSBarry Smith 
2200f1580f4eSBarry Smith    Should the functionality of `PCFieldSplitSchurGetSubKSP()` and `PCFieldSplitGetSubKSP()` be merged?
2201f1580f4eSBarry Smith 
2202f1580f4eSBarry Smith    The Fortran interface should be modernized to return directly the array of values.
2203f1580f4eSBarry Smith 
2204f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()`, `PCFieldSplitGetSubKSP()`
2205285fb4e2SStefano Zampini @*/
22069371c9d4SSatish Balay PetscErrorCode PCFieldSplitSchurGetSubKSP(PC pc, PetscInt *n, KSP *subksp[]) {
2207285fb4e2SStefano Zampini   PetscFunctionBegin;
2208285fb4e2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2209285fb4e2SStefano Zampini   if (n) PetscValidIntPointer(n, 2);
2210cac4c232SBarry Smith   PetscUseMethod(pc, "PCFieldSplitSchurGetSubKSP_C", (PC, PetscInt *, KSP **), (pc, n, subksp));
2211285fb4e2SStefano Zampini   PetscFunctionReturn(0);
2212285fb4e2SStefano Zampini }
2213285fb4e2SStefano Zampini 
2214e69d4d44SBarry Smith /*@
2215a1e3cbf9SBarry Smith     PCFieldSplitSetSchurPre -  Indicates from what operator the preconditioner is constructucted for the Schur complement.
2216a1e3cbf9SBarry Smith       The default is the A11 matrix.
2217e69d4d44SBarry Smith 
2218f1580f4eSBarry Smith     Collective on pc
2219e69d4d44SBarry Smith 
2220e69d4d44SBarry Smith     Input Parameters:
2221e69d4d44SBarry Smith +   pc      - the preconditioner context
2222f1580f4eSBarry Smith .   ptype   - which matrix to use for preconditioning the Schur complement: `PC_FIELDSPLIT_SCHUR_PRE_A11` (default),
2223f1580f4eSBarry Smith               `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PC_FIELDSPLIT_SCHUR_PRE_USER`,
2224f1580f4eSBarry Smith               `PC_FIELDSPLIT_SCHUR_PRE_SELFP`, and `PC_FIELDSPLIT_SCHUR_PRE_FULL`
22250298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
2226084e4875SJed Brown 
2227f1580f4eSBarry Smith     Options Database Keys:
2228a1e3cbf9SBarry Smith +    -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments
2229a1e3cbf9SBarry Smith -    -fieldsplit_1_pc_type <pctype> - the preconditioner algorithm that is used to construct the preconditioner from the operator
2230e69d4d44SBarry Smith 
2231fd1303e9SJungho Lee     Notes:
2232f1580f4eSBarry Smith     If ptype is
2233f1580f4eSBarry Smith +     a11 - the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
2234f1580f4eSBarry Smith      matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
2235f1580f4eSBarry Smith .     self - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
2236f1580f4eSBarry Smith           The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
2237f1580f4eSBarry Smith           preconditioner
2238f1580f4eSBarry Smith .     user - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
2239f1580f4eSBarry Smith           to this function).
2240f1580f4eSBarry Smith .     selfp - the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
2241f1580f4eSBarry Smith           This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
2242f1580f4eSBarry Smith           lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
2243f1580f4eSBarry Smith -     full - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation
2244f1580f4eSBarry Smith       computed internally by `PCFIELDSPLIT` (this is expensive)
2245f1580f4eSBarry Smith       useful mostly as a test that the Schur complement approach can work for your problem
2246fd1303e9SJungho Lee 
2247e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
2248fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
2249a7476a74SDmitry Karpeev     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
2250fd1303e9SJungho Lee 
2251e69d4d44SBarry Smith     Level: intermediate
2252e69d4d44SBarry Smith 
2253db781477SPatrick Sanan .seealso: `PCFieldSplitGetSchurPre()`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`,
2254f1580f4eSBarry Smith           `MatSchurComplementSetAinvType()`, `PCLSC`,
2255f1580f4eSBarry Smith           `PCFieldSplitSchurPreType`
2256e69d4d44SBarry Smith @*/
22579371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetSchurPre(PC pc, PCFieldSplitSchurPreType ptype, Mat pre) {
2258e69d4d44SBarry Smith   PetscFunctionBegin;
22590700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2260cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetSchurPre_C", (PC, PCFieldSplitSchurPreType, Mat), (pc, ptype, pre));
2261e69d4d44SBarry Smith   PetscFunctionReturn(0);
2262e69d4d44SBarry Smith }
2263686bed4dSStefano Zampini 
22649371c9d4SSatish Balay PetscErrorCode PCFieldSplitSchurPrecondition(PC pc, PCFieldSplitSchurPreType ptype, Mat pre) {
22659371c9d4SSatish Balay   return PCFieldSplitSetSchurPre(pc, ptype, pre);
22669371c9d4SSatish Balay } /* Deprecated name */
2267e69d4d44SBarry Smith 
226837a82bf0SJed Brown /*@
226937a82bf0SJed Brown     PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
2270f1580f4eSBarry Smith     preconditioned.  See `PCFieldSplitSetSchurPre()` for details.
227137a82bf0SJed Brown 
2272f1580f4eSBarry Smith     Logically Collective on pc
227337a82bf0SJed Brown 
2274f899ff85SJose E. Roman     Input Parameter:
227537a82bf0SJed Brown .   pc      - the preconditioner context
227637a82bf0SJed Brown 
227737a82bf0SJed Brown     Output Parameters:
2278f1580f4eSBarry Smith +   ptype   - which matrix to use for preconditioning the Schur complement: `PC_FIELDSPLIT_SCHUR_PRE_A11`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PC_FIELDSPLIT_PRE_USER`
2279f1580f4eSBarry Smith -   userpre - matrix to use for preconditioning (with `PC_FIELDSPLIT_PRE_USER`), or NULL
228037a82bf0SJed Brown 
228137a82bf0SJed Brown     Level: intermediate
228237a82bf0SJed Brown 
2283f1580f4eSBarry Smith .seealso: `PCFieldSplitSetSchurPre()`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`, `PCLSC`,
2284f1580f4eSBarry Smith           `PCFieldSplitSchurPreType`
228537a82bf0SJed Brown @*/
22869371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetSchurPre(PC pc, PCFieldSplitSchurPreType *ptype, Mat *pre) {
228737a82bf0SJed Brown   PetscFunctionBegin;
228837a82bf0SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2289cac4c232SBarry Smith   PetscUseMethod(pc, "PCFieldSplitGetSchurPre_C", (PC, PCFieldSplitSchurPreType *, Mat *), (pc, ptype, pre));
2290e69d4d44SBarry Smith   PetscFunctionReturn(0);
2291e69d4d44SBarry Smith }
2292e69d4d44SBarry Smith 
229345e7fc46SDmitry Karpeev /*@
2294f1580f4eSBarry Smith     PCFieldSplitSchurGetS -  extract the `MATSCHURCOMPLEMENT` object used by this `PCFIELDSPLIT` in case it needs to be configured separately
229545e7fc46SDmitry Karpeev 
229645e7fc46SDmitry Karpeev     Not collective
229745e7fc46SDmitry Karpeev 
229845e7fc46SDmitry Karpeev     Input Parameter:
229945e7fc46SDmitry Karpeev .   pc      - the preconditioner context
230045e7fc46SDmitry Karpeev 
230145e7fc46SDmitry Karpeev     Output Parameter:
2302470b340bSDmitry Karpeev .   S       - the Schur complement matrix
230345e7fc46SDmitry Karpeev 
2304f1580f4eSBarry Smith     Note:
2305f1580f4eSBarry Smith     This matrix should not be destroyed using `MatDestroy()`; rather, use `PCFieldSplitSchurRestoreS()`.
230645e7fc46SDmitry Karpeev 
230745e7fc46SDmitry Karpeev     Level: advanced
230845e7fc46SDmitry Karpeev 
2309f1580f4eSBarry Smith .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSchurPreType`, `PCFieldSplitSetSchurPre()`, `MATSCHURCOMPLEMENT`, `PCFieldSplitSchurRestoreS()`,
2310f1580f4eSBarry Smith           `MatCreateSchurComplement()`, `MatSchurComplementGetKSP()`, `MatSchurComplementComputeExplicitOperator()`, `MatGetSchurComplement()`
231145e7fc46SDmitry Karpeev @*/
23129371c9d4SSatish Balay PetscErrorCode PCFieldSplitSchurGetS(PC pc, Mat *S) {
231345e7fc46SDmitry Karpeev   const char    *t;
231445e7fc46SDmitry Karpeev   PetscBool      isfs;
231545e7fc46SDmitry Karpeev   PC_FieldSplit *jac;
231645e7fc46SDmitry Karpeev 
231745e7fc46SDmitry Karpeev   PetscFunctionBegin;
231845e7fc46SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
23199566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetType((PetscObject)pc, &t));
23209566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(t, PCFIELDSPLIT, &isfs));
232128b400f6SJacob Faibussowitsch   PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PC of type PCFIELDSPLIT, got %s instead", t);
232245e7fc46SDmitry Karpeev   jac = (PC_FieldSplit *)pc->data;
232363a3b9bcSJacob Faibussowitsch   PetscCheck(jac->type == PC_COMPOSITE_SCHUR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PCFIELDSPLIT of type SCHUR, got %d instead", jac->type);
2324470b340bSDmitry Karpeev   if (S) *S = jac->schur;
232545e7fc46SDmitry Karpeev   PetscFunctionReturn(0);
232645e7fc46SDmitry Karpeev }
232745e7fc46SDmitry Karpeev 
2328470b340bSDmitry Karpeev /*@
2329f1580f4eSBarry Smith     PCFieldSplitSchurRestoreS -  returns the `MATSCHURCOMPLEMENT` matrix used by this `PC`
2330470b340bSDmitry Karpeev 
2331470b340bSDmitry Karpeev     Not collective
2332470b340bSDmitry Karpeev 
2333470b340bSDmitry Karpeev     Input Parameters:
2334470b340bSDmitry Karpeev +   pc      - the preconditioner context
2335a2b725a8SWilliam Gropp -   S       - the Schur complement matrix
2336470b340bSDmitry Karpeev 
2337470b340bSDmitry Karpeev     Level: advanced
2338470b340bSDmitry Karpeev 
2339db781477SPatrick Sanan .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSchurPreType`, `PCFieldSplitSetSchurPre()`, `MatSchurComplement`, `PCFieldSplitSchurGetS()`
2340470b340bSDmitry Karpeev @*/
23419371c9d4SSatish Balay PetscErrorCode PCFieldSplitSchurRestoreS(PC pc, Mat *S) {
2342470b340bSDmitry Karpeev   const char    *t;
2343470b340bSDmitry Karpeev   PetscBool      isfs;
2344470b340bSDmitry Karpeev   PC_FieldSplit *jac;
2345470b340bSDmitry Karpeev 
2346470b340bSDmitry Karpeev   PetscFunctionBegin;
2347470b340bSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
23489566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetType((PetscObject)pc, &t));
23499566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(t, PCFIELDSPLIT, &isfs));
235028b400f6SJacob Faibussowitsch   PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PC of type PCFIELDSPLIT, got %s instead", t);
2351470b340bSDmitry Karpeev   jac = (PC_FieldSplit *)pc->data;
235263a3b9bcSJacob Faibussowitsch   PetscCheck(jac->type == PC_COMPOSITE_SCHUR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PCFIELDSPLIT of type SCHUR, got %d instead", jac->type);
23537827d75bSBarry Smith   PetscCheck(S && (*S == jac->schur), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatSchurComplement restored is not the same as gotten");
2354470b340bSDmitry Karpeev   PetscFunctionReturn(0);
2355470b340bSDmitry Karpeev }
2356470b340bSDmitry Karpeev 
23579371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetSchurPre_FieldSplit(PC pc, PCFieldSplitSchurPreType ptype, Mat pre) {
2358e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2359e69d4d44SBarry Smith 
2360e69d4d44SBarry Smith   PetscFunctionBegin;
2361084e4875SJed Brown   jac->schurpre = ptype;
2362a7476a74SDmitry Karpeev   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
23639566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->schur_user));
2364084e4875SJed Brown     jac->schur_user = pre;
23659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)jac->schur_user));
2366084e4875SJed Brown   }
2367e69d4d44SBarry Smith   PetscFunctionReturn(0);
2368e69d4d44SBarry Smith }
2369e69d4d44SBarry Smith 
23709371c9d4SSatish Balay static PetscErrorCode PCFieldSplitGetSchurPre_FieldSplit(PC pc, PCFieldSplitSchurPreType *ptype, Mat *pre) {
237137a82bf0SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
237237a82bf0SJed Brown 
237337a82bf0SJed Brown   PetscFunctionBegin;
237437a82bf0SJed Brown   *ptype = jac->schurpre;
237537a82bf0SJed Brown   *pre   = jac->schur_user;
237637a82bf0SJed Brown   PetscFunctionReturn(0);
237737a82bf0SJed Brown }
237837a82bf0SJed Brown 
2379ab1df9f5SJed Brown /*@
23800ffb0e17SBarry Smith     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain in the preconditioner
2381ab1df9f5SJed Brown 
2382f1580f4eSBarry Smith     Collective on pc
2383ab1df9f5SJed Brown 
2384ab1df9f5SJed Brown     Input Parameters:
2385ab1df9f5SJed Brown +   pc  - the preconditioner context
2386f1580f4eSBarry Smith -   ftype - which blocks of factorization to retain, `PC_FIELDSPLIT_SCHUR_FACT_FULL` is default
2387ab1df9f5SJed Brown 
2388f1580f4eSBarry Smith     Options Database Key:
2389147403d9SBarry Smith .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - default is full
2390ab1df9f5SJed Brown 
2391ab1df9f5SJed Brown     Level: intermediate
2392ab1df9f5SJed Brown 
2393ab1df9f5SJed Brown     Notes:
2394ab1df9f5SJed Brown     The FULL factorization is
2395ab1df9f5SJed Brown 
2396147403d9SBarry Smith .vb
2397147403d9SBarry Smith    (A   B)  = (1       0) (A   0) (1  Ainv*B)  = L D U
2398147403d9SBarry Smith    (C   E)    (C*Ainv  1) (0   S) (0     1)
2399147403d9SBarry Smith .vb
24000ffb0e17SBarry 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,
2401c096484dSStefano 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().
2402ab1df9f5SJed Brown 
2403147403d9SBarry Smith     If A and S are solved exactly
2404147403d9SBarry Smith .vb
2405147403d9SBarry Smith       *) FULL factorization is a direct solver.
2406f1580f4eSBarry 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.
2407f1580f4eSBarry 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.
2408147403d9SBarry Smith .ve
2409ab1df9f5SJed Brown 
2410f1580f4eSBarry Smith     If the iteration count is very low, consider using `KSPFGMRES` or `KSPGCR` which can use one less preconditioner
24110ffb0e17SBarry 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.
24120ffb0e17SBarry Smith 
2413f1580f4eSBarry Smith     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with `KSPMINRES`.
24140ffb0e17SBarry Smith 
2415f1580f4eSBarry Smith     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).
2416ab1df9f5SJed Brown 
2417ab1df9f5SJed Brown     References:
2418606c0280SSatish Balay +   * - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
2419606c0280SSatish Balay -   * - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).
2420ab1df9f5SJed Brown 
2421db781477SPatrick Sanan .seealso: `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`, `PCFieldSplitSetSchurScale()`
2422ab1df9f5SJed Brown @*/
24239371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetSchurFactType(PC pc, PCFieldSplitSchurFactType ftype) {
2424ab1df9f5SJed Brown   PetscFunctionBegin;
2425ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2426cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetSchurFactType_C", (PC, PCFieldSplitSchurFactType), (pc, ftype));
2427ab1df9f5SJed Brown   PetscFunctionReturn(0);
2428ab1df9f5SJed Brown }
2429ab1df9f5SJed Brown 
24309371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc, PCFieldSplitSchurFactType ftype) {
2431ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2432ab1df9f5SJed Brown 
2433ab1df9f5SJed Brown   PetscFunctionBegin;
2434ab1df9f5SJed Brown   jac->schurfactorization = ftype;
2435ab1df9f5SJed Brown   PetscFunctionReturn(0);
2436ab1df9f5SJed Brown }
2437ab1df9f5SJed Brown 
2438c096484dSStefano Zampini /*@
2439f1580f4eSBarry Smith     PCFieldSplitSetSchurScale -  Controls the sign flip of S for `PC_FIELDSPLIT_SCHUR_FACT_DIAG`.
2440c096484dSStefano Zampini 
2441f1580f4eSBarry Smith     Collective on pc
2442c096484dSStefano Zampini 
2443c096484dSStefano Zampini     Input Parameters:
2444c096484dSStefano Zampini +   pc    - the preconditioner context
2445c096484dSStefano Zampini -   scale - scaling factor for the Schur complement
2446c096484dSStefano Zampini 
2447f1580f4eSBarry Smith     Options Database Key:
2448c096484dSStefano Zampini .     -pc_fieldsplit_schur_scale - default is -1.0
2449c096484dSStefano Zampini 
2450c096484dSStefano Zampini     Level: intermediate
2451c096484dSStefano Zampini 
2452f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurFactType`, `PCFieldSplitSetSchurScale()`, `PCFieldSplitSetSchurFactType()`
2453c096484dSStefano Zampini @*/
24549371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetSchurScale(PC pc, PetscScalar scale) {
2455c096484dSStefano Zampini   PetscFunctionBegin;
2456c096484dSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2457c096484dSStefano Zampini   PetscValidLogicalCollectiveScalar(pc, scale, 2);
2458cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetSchurScale_C", (PC, PetscScalar), (pc, scale));
2459c096484dSStefano Zampini   PetscFunctionReturn(0);
2460c096484dSStefano Zampini }
2461c096484dSStefano Zampini 
24629371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc, PetscScalar scale) {
2463c096484dSStefano Zampini   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2464c096484dSStefano Zampini 
2465c096484dSStefano Zampini   PetscFunctionBegin;
2466c096484dSStefano Zampini   jac->schurscale = scale;
2467c096484dSStefano Zampini   PetscFunctionReturn(0);
2468c096484dSStefano Zampini }
2469c096484dSStefano Zampini 
247030ad9308SMatthew Knepley /*@C
24718c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
247230ad9308SMatthew Knepley 
2473f1580f4eSBarry Smith    Collective on pc
247430ad9308SMatthew Knepley 
247530ad9308SMatthew Knepley    Input Parameter:
247630ad9308SMatthew Knepley .  pc - the preconditioner context
247730ad9308SMatthew Knepley 
247830ad9308SMatthew Knepley    Output Parameters:
2479a04f6461SBarry Smith +  A00 - the (0,0) block
2480a04f6461SBarry Smith .  A01 - the (0,1) block
2481a04f6461SBarry Smith .  A10 - the (1,0) block
2482a04f6461SBarry Smith -  A11 - the (1,1) block
248330ad9308SMatthew Knepley 
248430ad9308SMatthew Knepley    Level: advanced
248530ad9308SMatthew Knepley 
2486f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `MatSchurComplementGetSubMatrices()`, `MatSchurComplementSetSubMatrices()`
248730ad9308SMatthew Knepley @*/
24889371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc, Mat *A00, Mat *A01, Mat *A10, Mat *A11) {
248930ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
249030ad9308SMatthew Knepley 
249130ad9308SMatthew Knepley   PetscFunctionBegin;
24920700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
249308401ef6SPierre Jolivet   PetscCheck(jac->type == PC_COMPOSITE_SCHUR, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2494a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
2495a04f6461SBarry Smith   if (A01) *A01 = jac->B;
2496a04f6461SBarry Smith   if (A10) *A10 = jac->C;
2497a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
249830ad9308SMatthew Knepley   PetscFunctionReturn(0);
249930ad9308SMatthew Knepley }
250030ad9308SMatthew Knepley 
2501b09e027eSCarola Kruse /*@
2502f1580f4eSBarry Smith     PCFieldSplitSetGKBTol -  Sets the solver tolerance for the generalized Golub-Kahan bidiagonalization preconditioner in `PCFIELDSPLIT`
2503b09e027eSCarola Kruse 
2504f1580f4eSBarry Smith     Collective on pc
2505e071a0a4SCarola Kruse 
2506b09e027eSCarola Kruse     Input Parameters:
2507b09e027eSCarola Kruse +   pc        - the preconditioner context
2508b09e027eSCarola Kruse -   tolerance - the solver tolerance
2509b09e027eSCarola Kruse 
2510f1580f4eSBarry Smith     Options Database Key:
2511b09e027eSCarola Kruse .   -pc_fieldsplit_gkb_tol - default is 1e-5
2512b09e027eSCarola Kruse 
2513b09e027eSCarola Kruse     Level: intermediate
2514b09e027eSCarola Kruse 
2515f1580f4eSBarry Smith     Note:
2516f1580f4eSBarry Smith     The generalized GKB algorithm uses a lower bound estimate of the error in energy norm as stopping criterion.
2517f1580f4eSBarry Smith     It stops once the lower bound estimate undershoots the required solver tolerance. Although the actual error might be bigger than
2518f1580f4eSBarry Smith     this estimate, the stopping criterion is satisfactory in practical cases [A13].
2519f1580f4eSBarry Smith 
2520f1580f4eSBarry Smith [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2521f1580f4eSBarry Smith 
2522db781477SPatrick Sanan .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetGKBDelay()`, `PCFieldSplitSetGKBNu()`, `PCFieldSplitSetGKBMaxit()`
2523b09e027eSCarola Kruse @*/
25249371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetGKBTol(PC pc, PetscReal tolerance) {
2525b09e027eSCarola Kruse   PetscFunctionBegin;
2526b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2527de482cd7SCarola Kruse   PetscValidLogicalCollectiveReal(pc, tolerance, 2);
2528cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetGKBTol_C", (PC, PetscReal), (pc, tolerance));
2529b09e027eSCarola Kruse   PetscFunctionReturn(0);
2530b09e027eSCarola Kruse }
2531b09e027eSCarola Kruse 
25329371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetGKBTol_FieldSplit(PC pc, PetscReal tolerance) {
2533b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2534b09e027eSCarola Kruse 
2535b09e027eSCarola Kruse   PetscFunctionBegin;
2536b09e027eSCarola Kruse   jac->gkbtol = tolerance;
2537b09e027eSCarola Kruse   PetscFunctionReturn(0);
2538b09e027eSCarola Kruse }
2539b09e027eSCarola Kruse 
2540b09e027eSCarola Kruse /*@
2541f1580f4eSBarry Smith     PCFieldSplitSetGKBMaxit -  Sets the maximum number of iterations for the generalized Golub-Kahan bidiagonalization preconditioner in `PCFIELDSPLIT`
2542b09e027eSCarola Kruse 
2543f1580f4eSBarry Smith     Collective on pc
2544b09e027eSCarola Kruse 
2545b09e027eSCarola Kruse     Input Parameters:
2546b09e027eSCarola Kruse +   pc     - the preconditioner context
2547b09e027eSCarola Kruse -   maxit  - the maximum number of iterations
2548b09e027eSCarola Kruse 
2549f1580f4eSBarry Smith     Options Database Key:
2550b09e027eSCarola Kruse .     -pc_fieldsplit_gkb_maxit - default is 100
2551b09e027eSCarola Kruse 
2552b09e027eSCarola Kruse     Level: intermediate
2553b09e027eSCarola Kruse 
2554db781477SPatrick Sanan .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetGKBDelay()`, `PCFieldSplitSetGKBTol()`, `PCFieldSplitSetGKBNu()`
2555b09e027eSCarola Kruse @*/
25569371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetGKBMaxit(PC pc, PetscInt maxit) {
2557b09e027eSCarola Kruse   PetscFunctionBegin;
2558b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2559de482cd7SCarola Kruse   PetscValidLogicalCollectiveInt(pc, maxit, 2);
2560cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetGKBMaxit_C", (PC, PetscInt), (pc, maxit));
2561b09e027eSCarola Kruse   PetscFunctionReturn(0);
2562b09e027eSCarola Kruse }
2563b09e027eSCarola Kruse 
25649371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetGKBMaxit_FieldSplit(PC pc, PetscInt maxit) {
2565b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2566b09e027eSCarola Kruse 
2567b09e027eSCarola Kruse   PetscFunctionBegin;
2568b09e027eSCarola Kruse   jac->gkbmaxit = maxit;
2569b09e027eSCarola Kruse   PetscFunctionReturn(0);
2570b09e027eSCarola Kruse }
2571b09e027eSCarola Kruse 
2572b09e027eSCarola Kruse /*@
2573f1580f4eSBarry Smith     PCFieldSplitSetGKBDelay -  Sets the delay in the lower bound error estimate in the generalized Golub-Kahan bidiagonalization in `PCFIELDSPLIT`
2574e071a0a4SCarola Kruse     preconditioner.
2575b09e027eSCarola Kruse 
2576f1580f4eSBarry Smith     Collective on pc
2577b09e027eSCarola Kruse 
2578b09e027eSCarola Kruse     Input Parameters:
2579b09e027eSCarola Kruse +   pc     - the preconditioner context
2580b09e027eSCarola Kruse -   delay  - the delay window in the lower bound estimate
2581b09e027eSCarola Kruse 
2582f1580f4eSBarry Smith     Options Database Key:
2583b09e027eSCarola Kruse .     -pc_fieldsplit_gkb_delay - default is 5
2584b09e027eSCarola Kruse 
2585b09e027eSCarola Kruse     Level: intermediate
2586b09e027eSCarola Kruse 
2587f1580f4eSBarry Smith     Note:
2588f1580f4eSBarry Smith     The algorithm uses a lower bound estimate of the error in energy norm as stopping criterion. The lower bound of the error ||u-u^k||_H
2589f1580f4eSBarry Smith     is expressed as a truncated sum. The error at iteration k can only be measured at iteration (k + delay), and thus the algorithm needs
2590f1580f4eSBarry Smith     at least (delay + 1) iterations to stop. For more details on the generalized Golub-Kahan bidiagonalization method and its lower bound stopping criterion, please refer to
2591f1580f4eSBarry Smith 
2592f1580f4eSBarry Smith [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2593f1580f4eSBarry Smith 
2594db781477SPatrick Sanan .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetGKBNu()`, `PCFieldSplitSetGKBTol()`, `PCFieldSplitSetGKBMaxit()`
2595b09e027eSCarola Kruse @*/
25969371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetGKBDelay(PC pc, PetscInt delay) {
2597b09e027eSCarola Kruse   PetscFunctionBegin;
2598b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2599b09e027eSCarola Kruse   PetscValidLogicalCollectiveInt(pc, delay, 2);
2600cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetGKBDelay_C", (PC, PetscInt), (pc, delay));
2601b09e027eSCarola Kruse   PetscFunctionReturn(0);
2602b09e027eSCarola Kruse }
2603b09e027eSCarola Kruse 
26049371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetGKBDelay_FieldSplit(PC pc, PetscInt delay) {
2605b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2606b09e027eSCarola Kruse 
2607b09e027eSCarola Kruse   PetscFunctionBegin;
2608b09e027eSCarola Kruse   jac->gkbdelay = delay;
2609b09e027eSCarola Kruse   PetscFunctionReturn(0);
2610b09e027eSCarola Kruse }
2611b09e027eSCarola Kruse 
2612b09e027eSCarola Kruse /*@
2613f1580f4eSBarry Smith     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
2614f1580f4eSBarry Smith     in `PCFIELDSPLIT`
2615b09e027eSCarola Kruse 
2616f1580f4eSBarry Smith     Collective on pc
2617f1580f4eSBarry Smith 
2618f1580f4eSBarry Smith     Input Parameters:
2619f1580f4eSBarry Smith +   pc     - the preconditioner context
2620f1580f4eSBarry Smith -   nu     - the shift parameter
2621f1580f4eSBarry Smith 
2622f1580f4eSBarry Smith     Options Database Keys:
2623f1580f4eSBarry Smith .     -pc_fieldsplit_gkb_nu - default is 1
2624f1580f4eSBarry Smith 
2625f1580f4eSBarry Smith     Level: intermediate
2626b09e027eSCarola Kruse 
2627b09e027eSCarola Kruse     Notes:
2628e071a0a4SCarola Kruse     This shift is in general done to obtain better convergence properties for the outer loop of the algorithm. This is often achieved by chosing nu sufficiently big. However,
2629b09e027eSCarola 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
2630b09e027eSCarola Kruse     necessary to find a good balance in between the convergence of the inner and outer loop.
2631b09e027eSCarola Kruse 
2632b09e027eSCarola 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.
2633b09e027eSCarola Kruse 
2634b09e027eSCarola Kruse [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2635b09e027eSCarola Kruse 
2636db781477SPatrick Sanan .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetGKBDelay()`, `PCFieldSplitSetGKBTol()`, `PCFieldSplitSetGKBMaxit()`
2637b09e027eSCarola Kruse @*/
26389371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetGKBNu(PC pc, PetscReal nu) {
2639b09e027eSCarola Kruse   PetscFunctionBegin;
2640b09e027eSCarola Kruse   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2641de482cd7SCarola Kruse   PetscValidLogicalCollectiveReal(pc, nu, 2);
2642cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetGKBNu_C", (PC, PetscReal), (pc, nu));
2643b09e027eSCarola Kruse   PetscFunctionReturn(0);
2644b09e027eSCarola Kruse }
2645b09e027eSCarola Kruse 
26469371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetGKBNu_FieldSplit(PC pc, PetscReal nu) {
2647b09e027eSCarola Kruse   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2648b09e027eSCarola Kruse 
2649b09e027eSCarola Kruse   PetscFunctionBegin;
2650b09e027eSCarola Kruse   jac->gkbnu = nu;
2651b09e027eSCarola Kruse   PetscFunctionReturn(0);
2652b09e027eSCarola Kruse }
2653b09e027eSCarola Kruse 
26549371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc, PCCompositeType type) {
265579416396SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
265679416396SBarry Smith 
265779416396SBarry Smith   PetscFunctionBegin;
265879416396SBarry Smith   jac->type = type;
26592e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", NULL));
26602e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurPre_C", NULL));
26612e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSchurPre_C", NULL));
26622e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurFactType_C", NULL));
26632e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurScale_C", NULL));
26642e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBTol_C", NULL));
26652e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBMaxit_C", NULL));
26662e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBNu_C", NULL));
26672e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBDelay_C", NULL));
2668e071a0a4SCarola Kruse 
26693b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
26703b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
26713b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
26722fa5cd67SKarl Rupp 
26739566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit_Schur));
26749566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurPre_C", PCFieldSplitSetSchurPre_FieldSplit));
26759566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSchurPre_C", PCFieldSplitGetSchurPre_FieldSplit));
26769566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurFactType_C", PCFieldSplitSetSchurFactType_FieldSplit));
26779566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurScale_C", PCFieldSplitSetSchurScale_FieldSplit));
2678a51937d4SCarola Kruse   } else if (type == PC_COMPOSITE_GKB) {
2679a51937d4SCarola Kruse     pc->ops->apply = PCApply_FieldSplit_GKB;
2680a51937d4SCarola Kruse     pc->ops->view  = PCView_FieldSplit_GKB;
2681e69d4d44SBarry Smith 
26829566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit));
26839566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBTol_C", PCFieldSplitSetGKBTol_FieldSplit));
26849566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBMaxit_C", PCFieldSplitSetGKBMaxit_FieldSplit));
26859566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBNu_C", PCFieldSplitSetGKBNu_FieldSplit));
26869566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBDelay_C", PCFieldSplitSetGKBDelay_FieldSplit));
26873b224e63SBarry Smith   } else {
26883b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
26893b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
26902fa5cd67SKarl Rupp 
26919566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit));
26923b224e63SBarry Smith   }
269379416396SBarry Smith   PetscFunctionReturn(0);
269479416396SBarry Smith }
269579416396SBarry Smith 
26969371c9d4SSatish Balay static PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc, PetscInt bs) {
269751f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
269851f519a2SBarry Smith 
269951f519a2SBarry Smith   PetscFunctionBegin;
270063a3b9bcSJacob Faibussowitsch   PetscCheck(bs >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %" PetscInt_FMT, bs);
27012472a847SBarry Smith   PetscCheck(jac->bs <= 0 || jac->bs == bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change fieldsplit blocksize from %" PetscInt_FMT " to %" PetscInt_FMT " after it has been set", jac->bs, bs);
270251f519a2SBarry Smith   jac->bs = bs;
270351f519a2SBarry Smith   PetscFunctionReturn(0);
270451f519a2SBarry Smith }
270551f519a2SBarry Smith 
27069371c9d4SSatish Balay static PetscErrorCode PCSetCoordinates_FieldSplit(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[]) {
27075ddf11f8SNicolas Barnafi   PC_FieldSplit    *jac           = (PC_FieldSplit *)pc->data;
27085ddf11f8SNicolas Barnafi   PC_FieldSplitLink ilink_current = jac->head;
27095ddf11f8SNicolas Barnafi   IS                is_owned;
27105ddf11f8SNicolas Barnafi 
27115ddf11f8SNicolas Barnafi   PetscFunctionBegin;
27125ddf11f8SNicolas Barnafi   jac->coordinates_set = PETSC_TRUE; // Internal flag
27139566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(pc->mat, &is_owned, PETSC_NULL));
27145ddf11f8SNicolas Barnafi 
27155ddf11f8SNicolas Barnafi   while (ilink_current) {
27165ddf11f8SNicolas Barnafi     // For each IS, embed it to get local coords indces
27175ddf11f8SNicolas Barnafi     IS              is_coords;
27185ddf11f8SNicolas Barnafi     PetscInt        ndofs_block;
27195ddf11f8SNicolas Barnafi     const PetscInt *block_dofs_enumeration; // Numbering of the dofs relevant to the current block
27205ddf11f8SNicolas Barnafi 
27215ddf11f8SNicolas Barnafi     // Setting drop to true for safety. It should make no difference.
27229566063dSJacob Faibussowitsch     PetscCall(ISEmbed(ilink_current->is, is_owned, PETSC_TRUE, &is_coords));
27239566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(is_coords, &ndofs_block));
27249566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is_coords, &block_dofs_enumeration));
27255ddf11f8SNicolas Barnafi 
27265ddf11f8SNicolas Barnafi     // Allocate coordinates vector and set it directly
27279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ndofs_block * dim, &(ilink_current->coords)));
27285ddf11f8SNicolas Barnafi     for (PetscInt dof = 0; dof < ndofs_block; ++dof) {
2729ad540459SPierre Jolivet       for (PetscInt d = 0; d < dim; ++d) (ilink_current->coords)[dim * dof + d] = coords[dim * block_dofs_enumeration[dof] + d];
27305ddf11f8SNicolas Barnafi     }
27315ddf11f8SNicolas Barnafi     ilink_current->dim   = dim;
27325ddf11f8SNicolas Barnafi     ilink_current->ndofs = ndofs_block;
27339566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is_coords, &block_dofs_enumeration));
27349566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_coords));
27355ddf11f8SNicolas Barnafi     ilink_current = ilink_current->next;
27365ddf11f8SNicolas Barnafi   }
27379566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_owned));
27385ddf11f8SNicolas Barnafi   PetscFunctionReturn(0);
27395ddf11f8SNicolas Barnafi }
27405ddf11f8SNicolas Barnafi 
2741bc08b0f1SBarry Smith /*@
2742f1580f4eSBarry Smith    PCFieldSplitSetType - Sets the type, `PCCompositeType`, of a `PCFIELDSPLIT`
274379416396SBarry Smith 
2744f1580f4eSBarry Smith    Collective on pc
274579416396SBarry Smith 
2746d8d19677SJose E. Roman    Input Parameters:
2747a2b725a8SWilliam Gropp +  pc - the preconditioner context
2748f1580f4eSBarry Smith -  type - `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE` (default), `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR`
274979416396SBarry Smith 
275079416396SBarry Smith    Options Database Key:
2751a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
275279416396SBarry Smith 
2753b02e2d75SMatthew G Knepley    Level: Intermediate
275479416396SBarry Smith 
2755f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCCompositeType`, `PCCompositeGetType()`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`,
2756f1580f4eSBarry Smith           `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR`
275779416396SBarry Smith @*/
27589371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetType(PC pc, PCCompositeType type) {
275979416396SBarry Smith   PetscFunctionBegin;
27600700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2761cac4c232SBarry Smith   PetscTryMethod(pc, "PCFieldSplitSetType_C", (PC, PCCompositeType), (pc, type));
276279416396SBarry Smith   PetscFunctionReturn(0);
276379416396SBarry Smith }
276479416396SBarry Smith 
2765b02e2d75SMatthew G Knepley /*@
2766f1580f4eSBarry Smith   PCFieldSplitGetType - Gets the type, `PCCompositeType`, of a `PCFIELDSPLIT`
2767b02e2d75SMatthew G Knepley 
2768b02e2d75SMatthew G Knepley   Not collective
2769b02e2d75SMatthew G Knepley 
2770b02e2d75SMatthew G Knepley   Input Parameter:
2771b02e2d75SMatthew G Knepley . pc - the preconditioner context
2772b02e2d75SMatthew G Knepley 
2773b02e2d75SMatthew G Knepley   Output Parameter:
2774f1580f4eSBarry Smith . type - `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE` (default), `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR`
2775b02e2d75SMatthew G Knepley 
2776b02e2d75SMatthew G Knepley   Level: Intermediate
2777b02e2d75SMatthew G Knepley 
2778f1580f4eSBarry Smith .seealso: `PCCompositeSetType()`, `PCFIELDSPLIT`, `PCCompositeType`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`,
2779f1580f4eSBarry Smith           `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR`
2780b02e2d75SMatthew G Knepley @*/
27819371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) {
2782b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
2783b02e2d75SMatthew G Knepley 
2784b02e2d75SMatthew G Knepley   PetscFunctionBegin;
2785b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2786b02e2d75SMatthew G Knepley   PetscValidIntPointer(type, 2);
2787b02e2d75SMatthew G Knepley   *type = jac->type;
2788b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
2789b02e2d75SMatthew G Knepley }
2790b02e2d75SMatthew G Knepley 
27914ab8060aSDmitry Karpeev /*@
2792f1580f4eSBarry Smith    PCFieldSplitSetDMSplits - Flags whether `DMCreateFieldDecomposition()` should be used to define the splits in a `PCFIELDSPLIT`, whenever possible.
27934ab8060aSDmitry Karpeev 
2794f1580f4eSBarry Smith    Logically Collective on pc
27954ab8060aSDmitry Karpeev 
27964ab8060aSDmitry Karpeev    Input Parameters:
27974ab8060aSDmitry Karpeev +  pc   - the preconditioner context
2798f1580f4eSBarry Smith -  flg  - boolean indicating whether to use field splits defined by the `DM`
27994ab8060aSDmitry Karpeev 
28004ab8060aSDmitry Karpeev    Options Database Key:
2801f1580f4eSBarry Smith .  -pc_fieldsplit_dm_splits <bool> - use the field splits defined by the `DM`
28024ab8060aSDmitry Karpeev 
28034ab8060aSDmitry Karpeev    Level: Intermediate
28044ab8060aSDmitry Karpeev 
2805f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitGetDMSplits()`, `PCFieldSplitSetFields()`, `PCFieldsplitSetIS()`
28064ab8060aSDmitry Karpeev @*/
28079371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetDMSplits(PC pc, PetscBool flg) {
28084ab8060aSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
28094ab8060aSDmitry Karpeev   PetscBool      isfs;
28104ab8060aSDmitry Karpeev 
28114ab8060aSDmitry Karpeev   PetscFunctionBegin;
28124ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
28134ab8060aSDmitry Karpeev   PetscValidLogicalCollectiveBool(pc, flg, 2);
28149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs));
2815ad540459SPierre Jolivet   if (isfs) jac->dm_splits = flg;
28164ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
28174ab8060aSDmitry Karpeev }
28184ab8060aSDmitry Karpeev 
28194ab8060aSDmitry Karpeev /*@
2820f1580f4eSBarry Smith    PCFieldSplitGetDMSplits - Returns flag indicating whether `DMCreateFieldDecomposition()` should be used to define the splits in a `PCFIELDSPLIT`, whenever possible.
28214ab8060aSDmitry Karpeev 
28224ab8060aSDmitry Karpeev    Logically Collective
28234ab8060aSDmitry Karpeev 
28244ab8060aSDmitry Karpeev    Input Parameter:
28254ab8060aSDmitry Karpeev .  pc   - the preconditioner context
28264ab8060aSDmitry Karpeev 
28274ab8060aSDmitry Karpeev    Output Parameter:
2828f1580f4eSBarry Smith .  flg  - boolean indicating whether to use field splits defined by the `DM`
28294ab8060aSDmitry Karpeev 
28304ab8060aSDmitry Karpeev    Level: Intermediate
28314ab8060aSDmitry Karpeev 
2832f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetDMSplits()`, `PCFieldSplitSetFields()`, `PCFieldsplitSetIS()`
28334ab8060aSDmitry Karpeev @*/
28349371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetDMSplits(PC pc, PetscBool *flg) {
28354ab8060aSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
28364ab8060aSDmitry Karpeev   PetscBool      isfs;
28374ab8060aSDmitry Karpeev 
28384ab8060aSDmitry Karpeev   PetscFunctionBegin;
28394ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2840534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
28419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs));
28424ab8060aSDmitry Karpeev   if (isfs) {
28434ab8060aSDmitry Karpeev     if (flg) *flg = jac->dm_splits;
28444ab8060aSDmitry Karpeev   }
28454ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
28464ab8060aSDmitry Karpeev }
28474ab8060aSDmitry Karpeev 
28487b752e3dSPatrick Sanan /*@
2849f1580f4eSBarry Smith    PCFieldSplitGetDetectSaddlePoint - Returns flag indicating whether `PCFIELDSPLIT` will attempt to automatically determine fields based on zero diagonal entries.
28507b752e3dSPatrick Sanan 
28517b752e3dSPatrick Sanan    Logically Collective
28527b752e3dSPatrick Sanan 
28537b752e3dSPatrick Sanan    Input Parameter:
28547b752e3dSPatrick Sanan .  pc   - the preconditioner context
28557b752e3dSPatrick Sanan 
28567b752e3dSPatrick Sanan    Output Parameter:
28577b752e3dSPatrick Sanan .  flg  - boolean indicating whether to detect fields or not
28587b752e3dSPatrick Sanan 
28597b752e3dSPatrick Sanan    Level: Intermediate
28607b752e3dSPatrick Sanan 
2861db781477SPatrick Sanan .seealso: `PCFIELDSPLIT`, `PCFieldSplitSetDetectSaddlePoint()`
28627b752e3dSPatrick Sanan @*/
28639371c9d4SSatish Balay PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC pc, PetscBool *flg) {
28647b752e3dSPatrick Sanan   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
28657b752e3dSPatrick Sanan 
28667b752e3dSPatrick Sanan   PetscFunctionBegin;
28677b752e3dSPatrick Sanan   *flg = jac->detect;
28687b752e3dSPatrick Sanan   PetscFunctionReturn(0);
28697b752e3dSPatrick Sanan }
28707b752e3dSPatrick Sanan 
28717b752e3dSPatrick Sanan /*@
2872f1580f4eSBarry Smith    PCFieldSplitSetDetectSaddlePoint - Sets flag indicating whether `PCFIELDSPLIT` will attempt to automatically determine fields based on zero diagonal entries.
28737b752e3dSPatrick Sanan 
28747b752e3dSPatrick Sanan    Logically Collective
28757b752e3dSPatrick Sanan 
28767b752e3dSPatrick Sanan    Input Parameter:
28777b752e3dSPatrick Sanan .  pc   - the preconditioner context
28787b752e3dSPatrick Sanan 
28797b752e3dSPatrick Sanan    Output Parameter:
28807b752e3dSPatrick Sanan .  flg  - boolean indicating whether to detect fields or not
28817b752e3dSPatrick Sanan 
28827b752e3dSPatrick Sanan    Options Database Key:
2883147403d9SBarry Smith .  -pc_fieldsplit_detect_saddle_point <bool> - detect and use the saddle point
2884147403d9SBarry Smith 
2885f1580f4eSBarry Smith    Note:
2886f1580f4eSBarry Smith    Also sets the split type to `PC_COMPOSITE_SCHUR` (see `PCFieldSplitSetType()`) and the Schur preconditioner type to `PC_FIELDSPLIT_SCHUR_PRE_SELF` (see `PCFieldSplitSetSchurPre()`).
28877b752e3dSPatrick Sanan 
28887b752e3dSPatrick Sanan    Level: Intermediate
28897b752e3dSPatrick Sanan 
2890f1580f4eSBarry Smith .seealso: `PCFIELDSPLIT`, `PCFieldSplitGetDetectSaddlePoint()`, `PCFieldSplitSetType()`, `PCFieldSplitSetSchurPre()`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`
28917b752e3dSPatrick Sanan @*/
28929371c9d4SSatish Balay PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc, PetscBool flg) {
28937b752e3dSPatrick Sanan   PC_FieldSplit *jac = (PC_FieldSplit *)pc->data;
28947b752e3dSPatrick Sanan 
28957b752e3dSPatrick Sanan   PetscFunctionBegin;
28967b752e3dSPatrick Sanan   jac->detect = flg;
28977b752e3dSPatrick Sanan   if (jac->detect) {
28989566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR));
28999566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_SELF, NULL));
29007b752e3dSPatrick Sanan   }
29017b752e3dSPatrick Sanan   PetscFunctionReturn(0);
29027b752e3dSPatrick Sanan }
29037b752e3dSPatrick Sanan 
29040971522cSBarry Smith /*MC
2905a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2906f1580f4eSBarry Smith    collections of variables (that may overlap) called splits. See [the users manual section on "Solving Block Matrices"](sec_block_matrices) for more details.
29070971522cSBarry Smith 
2908f1580f4eSBarry Smith      To set options on the solvers for each block append `-fieldsplit_` to all the `PC`
2909de37d9f1SPatrick Sanan         options database keys. For example, `-fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1`
29100971522cSBarry Smith 
2911de37d9f1SPatrick Sanan      To set the options on the solvers separate for each block call `PCFieldSplitGetSubKSP()`
2912de37d9f1SPatrick Sanan          and set the options directly on the resulting `KSP` object
29130971522cSBarry Smith 
29140971522cSBarry Smith    Level: intermediate
29150971522cSBarry Smith 
291679416396SBarry Smith    Options Database Keys:
2917de37d9f1SPatrick Sanan +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the `%d`'th split
291881540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
2919de37d9f1SPatrick Sanan                               been supplied explicitly by `-pc_fieldsplit_%d_fields`
292081540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
2921a51937d4SCarola Kruse .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur,gkb> - type of relaxation or factorization splitting
2922de37d9f1SPatrick Sanan .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see `PCFieldSplitSetSchurPre()`
2923de37d9f1SPatrick Sanan .   -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - set factorization type when using `-pc_fieldsplit_type schur`; see `PCFieldSplitSetSchurFactType()`
2924fb6809a2SPatrick Sanan -   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero diagonal and uses Schur complement with no preconditioner as the solver
292579416396SBarry Smith 
2926de37d9f1SPatrick Sanan      Options prefixes for inner solvers when using the Schur complement preconditioner are `-fieldsplit_0_` and `-fieldsplit_1_` .
2927de37d9f1SPatrick Sanan      For all other solvers they are `-fieldsplit_%d_` for the `d`th field; use `-fieldsplit_` for all fields.
2928de37d9f1SPatrick Sanan      The options prefix for the inner solver when using the Golub-Kahan biadiagonalization preconditioner is `-fieldsplit_0_`
29295d4c12cdSJungho Lee 
2930c8a0d604SMatthew G Knepley    Notes:
2931f1580f4eSBarry Smith     Use `PCFieldSplitSetFields()` to set splits defined by "strided" entries and `PCFieldSplitSetIS()`
2932f1580f4eSBarry Smith      to define a split by an arbitrary collection of entries.
2933d32f9abdSBarry Smith 
2934f1580f4eSBarry Smith       If no splits are set the default is used. The splits are defined by entries strided by bs,
2935de37d9f1SPatrick Sanan       beginning at 0 then 1, etc to bs-1. The block size can be set with `PCFieldSplitSetBlockSize()`,
2936d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
2937de37d9f1SPatrick Sanan       to `KSPSetOperators()`/`PCSetOperators()`.
2938d32f9abdSBarry Smith 
2939de37d9f1SPatrick Sanan       For the Schur complement preconditioner if
29401d71051fSDmitry Karpeev 
2941de37d9f1SPatrick Sanan       ```{math}
2942de37d9f1SPatrick Sanan       J = \left[\begin{array}{cc} A_{00} & A_{01} \\ A_{10} & A_{11} \end{array}\right]
2943de37d9f1SPatrick Sanan       ```
2944e69d4d44SBarry Smith 
2945de37d9f1SPatrick Sanan       the preconditioner using `full` factorization is logically
2946de37d9f1SPatrick Sanan       ```{math}
2947de37d9f1SPatrick Sanan       \left[\begin{array}{cc} I & -\text{ksp}(A_{00}) \\ 0 & I \end{array}\right] \left[\begin{array}{cc} \text{inv}(A_{00}) & 0 \\ 0 & \text{ksp}(S) \end{array}\right] \left[\begin{array}{cc} I & 0 \\ -A_{10} \text{ksp}(A_{00}) & I \end{array}\right]
2948de37d9f1SPatrick Sanan       ```
2949de37d9f1SPatrick Sanan      where the action of $\text{inv}(A_{00})$ is applied using the KSP solver with prefix `-fieldsplit_0_`.  $S$ is the Schur complement
2950de37d9f1SPatrick Sanan      ```{math}
2951223e5b4fSPatrick Sanan      S = A_{11} - A_{10} \text{ksp}(A_{00}) A_{01}
2952de37d9f1SPatrick Sanan      ```
2953de37d9f1SPatrick Sanan      which is usually dense and not stored explicitly.  The action of $\text{ksp}(S)$ is computed using the KSP solver with prefix `-fieldsplit_splitname_` (where `splitname` was given
2954de37d9f1SPatrick Sanan      in providing the SECOND split or 1 if not given). For `PCFieldSplitGetSub\text{ksp}()` when field number is 0,
2955de37d9f1SPatrick Sanan      it returns the KSP associated with `-fieldsplit_0_` while field number 1 gives `-fieldsplit_1_` KSP. By default
2956de37d9f1SPatrick Sanan      $A_{11}$ is used to construct a preconditioner for $S$, use `PCFieldSplitSetSchurPre()` for all the possible ways to construct the preconditioner for $S$.
2957de37d9f1SPatrick Sanan 
2958de37d9f1SPatrick Sanan      The factorization type is set using `-pc_fieldsplit_schur_fact_type <diag, lower, upper, full>`. `full` is shown above,
2959de37d9f1SPatrick Sanan      `diag` gives
2960de37d9f1SPatrick Sanan       ```{math}
2961de37d9f1SPatrick Sanan       \left[\begin{array}{cc} \text{inv}(A_{00}) & 0 \\  0 & -\text{ksp}(S) \end{array}\right]
2962de37d9f1SPatrick Sanan       ```
2963de37d9f1SPatrick Sanan      Note that, slightly counter intuitively, there is a negative in front of the $\text{ksp}(S)$  so that the preconditioner is positive definite. For SPD matrices $J$, the sign flip
2964de37d9f1SPatrick Sanan      can be turned off with `PCFieldSplitSetSchurScale()` or by command line `-pc_fieldsplit_schur_scale 1.0`. The `lower` factorization is the inverse of
2965de37d9f1SPatrick Sanan       ```{math}
2966de37d9f1SPatrick Sanan       \left[\begin{array}{cc} A_{00} & 0 \\  A_{10} & S \end{array}\right]
2967de37d9f1SPatrick Sanan       ```
2968de37d9f1SPatrick Sanan      where the inverses of A_{00} and S are applied using KSPs. The upper factorization is the inverse of
2969de37d9f1SPatrick Sanan       ```{math}
2970de37d9f1SPatrick Sanan       \left[\begin{array}{cc} A_{00} & A_{01} \\  0 & S \end{array}\right]
2971de37d9f1SPatrick Sanan       ```
2972de37d9f1SPatrick Sanan      where again the inverses of $A_{00}$ and $S$ are applied using `KSP`s.
2973de37d9f1SPatrick Sanan 
2974de37d9f1SPatrick Sanan      If only one set of indices (one `IS`) is provided with `PCFieldSplitSetIS()` then the complement of that `IS`
2975edf189efSBarry Smith      is used automatically for a second block.
2976edf189efSBarry Smith 
2977f1580f4eSBarry Smith      The fieldsplit preconditioner cannot currently be used with the `MATBAIJ` or `MATSBAIJ` data formats if the blocksize is larger than 1.
2978f1580f4eSBarry Smith      Generally it should be used with the `MATAIJ` format.
2979ff218e97SBarry Smith 
2980ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
2981f1580f4eSBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling {cite}`Wesseling2009`.
2982f1580f4eSBarry Smith      One can also use `PCFIELDSPLIT`
2983ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
29840716a85fSBarry Smith 
2985de37d9f1SPatrick Sanan      See "A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations" {cite}`elman2008tcp`.
2986a6a584a2SBarry Smith 
2987de37d9f1SPatrick Sanan      The Constrained Pressure Preconditioner (CPR) can be implemented using `PCCOMPOSITE` with `PCGALERKIN`. CPR first solves an $R A P$ subsystem, updates the
2988de37d9f1SPatrick Sanan      residual on all variables (`PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE)`), and then applies a simple ILU like preconditioner on all the variables.
2989a51937d4SCarola Kruse 
2990de37d9f1SPatrick Sanan      The generalized Golub-Kahan bidiagonalization preconditioner (GKB) can be applied to symmetric $2 \times 2$ block matrices of the shape
2991de37d9f1SPatrick Sanan      ```{math}
2992de37d9f1SPatrick Sanan      \left[\begin{array}{cc} A_{00} & A_{01} \\ A_{01}' & 0 \end{array}\right]
2993de37d9f1SPatrick Sanan      ```
2994de37d9f1SPatrick Sanan      with $A_{00}$ positive semi-definite. The implementation follows {cite}`Arioli2013`. Therein, we choose $N := 1/\nu * I$ and the $(1,1)$-block of the matrix is modified to $H = _{A00} + \nu*A_{01}*A_{01}'$.
2995de37d9f1SPatrick Sanan      A linear system $Hx = b$ has to be solved in each iteration of the GKB algorithm. This solver is chosen with the option prefix `-fieldsplit_0_`.
2996a51937d4SCarola Kruse 
2997f1580f4eSBarry Smith      References:
2998f1580f4eSBarry Smith 
2999de37d9f1SPatrick Sanan      ```{bibliography}
3000de37d9f1SPatrick Sanan      :filter: docname in docnames
3001de37d9f1SPatrick Sanan      ```
3002de37d9f1SPatrick Sanan 
3003f1580f4eSBarry Smith    Developer Note:
3004f1580f4eSBarry Smith    The Schur complement functionality of `PCFIELDSPLIT` should likely be factored into its own `PC` thus simplify the implementation of the preconditioners and their
3005f1580f4eSBarry Smith    user API.
3006f1580f4eSBarry Smith 
3007db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCLSC`,
3008db781477SPatrick Sanan           `PCFieldSplitGetSubKSP()`, `PCFieldSplitSchurGetSubKSP()`, `PCFieldSplitSetFields()`,
3009db781477SPatrick Sanan           `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`, `PCFieldSplitSetSchurFactType()`,
3010db781477SPatrick Sanan           `MatSchurComplementSetAinvType()`, `PCFieldSplitSetSchurScale()`, `PCFieldSplitSetDetectSaddlePoint()`
30110971522cSBarry Smith M*/
30120971522cSBarry Smith 
30139371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc) {
30140971522cSBarry Smith   PC_FieldSplit *jac;
30150971522cSBarry Smith 
30160971522cSBarry Smith   PetscFunctionBegin;
3017*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
30182fa5cd67SKarl Rupp 
30190971522cSBarry Smith   jac->bs                 = -1;
30200971522cSBarry Smith   jac->nsplits            = 0;
30213e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
3022e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
3023c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
3024c096484dSStefano Zampini   jac->schurscale         = -1.0;
3025fbe7908bSJed Brown   jac->dm_splits          = PETSC_TRUE;
30267b752e3dSPatrick Sanan   jac->detect             = PETSC_FALSE;
3027a51937d4SCarola Kruse   jac->gkbtol             = 1e-5;
3028a51937d4SCarola Kruse   jac->gkbdelay           = 5;
3029a51937d4SCarola Kruse   jac->gkbnu              = 1;
3030a51937d4SCarola Kruse   jac->gkbmaxit           = 100;
3031de482cd7SCarola Kruse   jac->gkbmonitor         = PETSC_FALSE;
30325ddf11f8SNicolas Barnafi   jac->coordinates_set    = PETSC_FALSE;
303351f519a2SBarry Smith 
30340971522cSBarry Smith   pc->data = (void *)jac;
30350971522cSBarry Smith 
30360971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
3037421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
30380971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
3039574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
30400971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
30410971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
30420971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
30430a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
30440971522cSBarry Smith 
30459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSchurGetSubKSP_C", PCFieldSplitSchurGetSubKSP_FieldSplit));
30469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit));
30479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetFields_C", PCFieldSplitSetFields_FieldSplit));
30489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetIS_C", PCFieldSplitSetIS_FieldSplit));
30499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetType_C", PCFieldSplitSetType_FieldSplit));
30509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetBlockSize_C", PCFieldSplitSetBlockSize_FieldSplit));
30519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitRestrictIS_C", PCFieldSplitRestrictIS_FieldSplit));
30529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_FieldSplit));
30530971522cSBarry Smith   PetscFunctionReturn(0);
30540971522cSBarry Smith }
3055