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> 4f7cbcdf3SPierre Jolivet #include <petscdevice.h> 5f7cbcdf3SPierre Jolivet #if PetscDefined(HAVE_CUDA) 6f7cbcdf3SPierre Jolivet #include <petscdevice_cuda.h> 7f7cbcdf3SPierre Jolivet #endif 8f7cbcdf3SPierre Jolivet #if PetscDefined(HAVE_HIP) 9f7cbcdf3SPierre Jolivet #include <petscdevice_hip.h> 10f7cbcdf3SPierre Jolivet #endif 110971522cSBarry Smith 120a545947SLisandro Dalcin const char *const PCFieldSplitSchurPreTypes[] = {"SELF", "SELFP", "A11", "USER", "FULL", "PCFieldSplitSchurPreType", "PC_FIELDSPLIT_SCHUR_PRE_", NULL}; 130a545947SLisandro Dalcin const char *const PCFieldSplitSchurFactTypes[] = {"DIAG", "LOWER", "UPPER", "FULL", "PCFieldSplitSchurFactType", "PC_FIELDSPLIT_SCHUR_FACT_", NULL}; 14c5d2311dSJed Brown 159df09d43SBarry 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; 169df09d43SBarry Smith 170971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 180971522cSBarry Smith struct _PC_FieldSplitLink { 1969a612a9SBarry Smith KSP ksp; 20443836d0SMatthew G Knepley Vec x, y, z; 21db4c96c1SJed Brown char *splitname; 220971522cSBarry Smith PetscInt nfields; 235d4c12cdSJungho Lee PetscInt *fields, *fields_col; 241b9fc7fcSBarry Smith VecScatter sctx; 25f5f0d762SBarry Smith IS is, is_col; 2651f519a2SBarry Smith PC_FieldSplitLink next, previous; 279df09d43SBarry Smith PetscLogEvent event; 285ddf11f8SNicolas Barnafi 295ddf11f8SNicolas Barnafi /* Used only when setting coordinates with PCSetCoordinates */ 305ddf11f8SNicolas Barnafi PetscInt dim; 315ddf11f8SNicolas Barnafi PetscInt ndofs; 325ddf11f8SNicolas Barnafi PetscReal *coords; 330971522cSBarry Smith }; 340971522cSBarry Smith 350971522cSBarry Smith typedef struct { 3668dd23aaSBarry Smith PCCompositeType type; 37ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 38ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 3930ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 4030ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 4179416396SBarry Smith Vec *x, *y, w1, w2; 42519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 43519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 4430ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 45ace3abfcSBarry Smith PetscBool issetup; 462fa5cd67SKarl Rupp 4730ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 4830ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 4930ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 50443836d0SMatthew G Knepley Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */ 517addb90fSBarry Smith Mat schurp; /* Assembled approximation to S built by MatSchurComplement to be used as a matrix for constructing the preconditioner when solving with S */ 527addb90fSBarry Smith Mat schur_user; /* User-provided matrix for constructing the preconditioner for the Schur complement */ 537addb90fSBarry Smith PCFieldSplitSchurPreType schurpre; /* Determines which matrix is used for the Schur complement */ 54c9c6ffaaSJed Brown PCFieldSplitSchurFactType schurfactorization; 5530ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 56443836d0SMatthew G Knepley KSP kspupper; /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */ 57c096484dSStefano Zampini PetscScalar schurscale; /* Scaling factor for the Schur complement solution with DIAG factorization */ 58c096484dSStefano Zampini 59a51937d4SCarola Kruse /* Only used when Golub-Kahan bidiagonalization preconditioning is used */ 60a51937d4SCarola Kruse Mat H; /* The modified matrix H = A00 + nu*A01*A01' */ 61a51937d4SCarola Kruse PetscReal gkbtol; /* Stopping tolerance for lower bound estimate */ 62a51937d4SCarola Kruse PetscInt gkbdelay; /* The delay window for the stopping criterion */ 63a51937d4SCarola Kruse PetscReal gkbnu; /* Parameter for augmented Lagrangian H = A + nu*A01*A01' */ 64a51937d4SCarola Kruse PetscInt gkbmaxit; /* Maximum number of iterations for outer loop */ 65de482cd7SCarola Kruse PetscBool gkbmonitor; /* Monitor for gkb iterations and the lower bound error */ 66de482cd7SCarola Kruse PetscViewer gkbviewer; /* Viewer context for gkbmonitor */ 67e071a0a4SCarola Kruse Vec u, v, d, Hu; /* Work vectors for the GKB algorithm */ 68e071a0a4SCarola Kruse PetscScalar *vecz; /* Contains intermediate values, eg for lower bound */ 69a51937d4SCarola Kruse 7097bbdb24SBarry Smith PC_FieldSplitLink head; 716dbb499eSCian Wilson PetscBool isrestrict; /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */ 72c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 734ab8060aSDmitry Karpeev PetscBool dm_splits; /* Whether to use DMCreateFieldDecomposition() whenever possible */ 74c84da90fSDmitry Karpeev PetscBool diag_use_amat; /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */ 75c84da90fSDmitry Karpeev PetscBool offdiag_use_amat; /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */ 767b752e3dSPatrick Sanan PetscBool detect; /* Whether to form 2-way split by finding zero diagonal entries */ 775ddf11f8SNicolas Barnafi PetscBool coordinates_set; /* Whether PCSetCoordinates has been called */ 780971522cSBarry Smith } PC_FieldSplit; 790971522cSBarry Smith 8016913363SBarry Smith /* 81f1580f4eSBarry Smith Note: 8295452b02SPatrick Sanan there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 8316913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 8416913363SBarry Smith PC you could change this. 8516913363SBarry Smith */ 86084e4875SJed Brown 877addb90fSBarry Smith /* This helper is so that setting a user-provided matrix is orthogonal to choosing to use it. This way the 88084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 89d71ae5a4SJacob Faibussowitsch static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 90d71ae5a4SJacob Faibussowitsch { 91084e4875SJed Brown switch (jac->schurpre) { 92d71ae5a4SJacob Faibussowitsch case PC_FIELDSPLIT_SCHUR_PRE_SELF: 93d71ae5a4SJacob Faibussowitsch return jac->schur; 94d71ae5a4SJacob Faibussowitsch case PC_FIELDSPLIT_SCHUR_PRE_SELFP: 95d71ae5a4SJacob Faibussowitsch return jac->schurp; 96d71ae5a4SJacob Faibussowitsch case PC_FIELDSPLIT_SCHUR_PRE_A11: 97d71ae5a4SJacob Faibussowitsch return jac->pmat[1]; 98e74569cdSMatthew G. Knepley case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */ 99084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 100d71ae5a4SJacob Faibussowitsch default: 101d71ae5a4SJacob Faibussowitsch return jac->schur_user ? jac->schur_user : jac->pmat[1]; 102084e4875SJed Brown } 103084e4875SJed Brown } 104084e4875SJed Brown 1059804daf3SBarry Smith #include <petscdraw.h> 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_FieldSplit(PC pc, PetscViewer viewer) 107d71ae5a4SJacob Faibussowitsch { 1080971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 1099f196a02SMartin Diehl PetscBool isascii, isdraw; 1100971522cSBarry Smith PetscInt i, j; 1115a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 1120971522cSBarry Smith 1130971522cSBarry Smith PetscFunctionBegin; 1149f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1169f196a02SMartin Diehl if (isascii) { 1172c9966d7SBarry Smith if (jac->bs > 0) { 11863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FieldSplit with %s composition: total splits = %" PetscInt_FMT ", blocksize = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits, jac->bs)); 1192c9966d7SBarry Smith } else { 12063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FieldSplit with %s composition: total splits = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits)); 1212c9966d7SBarry Smith } 12248a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for blocks\n")); 12348a46eb9SPierre Jolivet if (jac->diag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for diagonal blocks\n")); 12448a46eb9SPierre Jolivet if (jac->offdiag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for off-diagonal blocks\n")); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Solver info for each split is in the following KSP objects:\n")); 1260971522cSBarry Smith for (i = 0; i < jac->nsplits; i++) { 1271ab39975SBarry Smith if (ilink->fields) { 12863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Fields ", i)); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 1305a9f2f41SSatish Balay for (j = 0; j < ilink->nfields; j++) { 13148a46eb9SPierre Jolivet if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ",")); 13263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, ilink->fields[j])); 1330971522cSBarry Smith } 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 1361ab39975SBarry Smith } else { 13763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Defined by IS\n", i)); 1381ab39975SBarry Smith } 1399566063dSJacob Faibussowitsch PetscCall(KSPView(ilink->ksp, viewer)); 1405a9f2f41SSatish Balay ilink = ilink->next; 1410971522cSBarry Smith } 1422fa5cd67SKarl Rupp } 1432fa5cd67SKarl Rupp 1442fa5cd67SKarl Rupp if (isdraw) { 145d9884438SBarry Smith PetscDraw draw; 146d9884438SBarry Smith PetscReal x, y, w, wd; 147d9884438SBarry Smith 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 1499566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 150d9884438SBarry Smith w = 2 * PetscMin(1.0 - x, x); 151d9884438SBarry Smith wd = w / (jac->nsplits + 1); 152d9884438SBarry Smith x = x - wd * (jac->nsplits - 1) / 2.0; 153d9884438SBarry Smith for (i = 0; i < jac->nsplits; i++) { 1549566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, y)); 1559566063dSJacob Faibussowitsch PetscCall(KSPView(ilink->ksp, viewer)); 1569566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 157d9884438SBarry Smith x += wd; 158d9884438SBarry Smith ilink = ilink->next; 159d9884438SBarry Smith } 1600971522cSBarry Smith } 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1620971522cSBarry Smith } 1630971522cSBarry Smith 164d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_FieldSplit_Schur(PC pc, PetscViewer viewer) 165d71ae5a4SJacob Faibussowitsch { 1663b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 1679f196a02SMartin Diehl PetscBool isascii, isdraw; 1683b224e63SBarry Smith PetscInt i, j; 1693b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 170a9908d51SBarry Smith MatSchurComplementAinvType atype; 1713b224e63SBarry Smith 1723b224e63SBarry Smith PetscFunctionBegin; 1739f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 1749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1759f196a02SMartin Diehl if (isascii) { 1762c9966d7SBarry Smith if (jac->bs > 0) { 17763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FieldSplit with Schur preconditioner, blocksize = %" PetscInt_FMT ", factorization %s\n", jac->bs, PCFieldSplitSchurFactTypes[jac->schurfactorization])); 1782c9966d7SBarry Smith } else { 1799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FieldSplit with Schur preconditioner, factorization %s\n", PCFieldSplitSchurFactTypes[jac->schurfactorization])); 1802c9966d7SBarry Smith } 18148a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for blocks\n")); 1823e8b8b31SMatthew G Knepley switch (jac->schurpre) { 183d71ae5a4SJacob Faibussowitsch case PC_FIELDSPLIT_SCHUR_PRE_SELF: 184d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Preconditioner for the Schur complement formed from S itself\n")); 185d71ae5a4SJacob Faibussowitsch break; 186a7476a74SDmitry Karpeev case PC_FIELDSPLIT_SCHUR_PRE_SELFP: 1877cf5f706SPierre Jolivet if (jac->schur) { 1889566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetAinvType(jac->schur, &atype)); 1899371c9d4SSatish 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 ")))); 1907cf5f706SPierre Jolivet } 191a9908d51SBarry Smith break; 192d71ae5a4SJacob Faibussowitsch case PC_FIELDSPLIT_SCHUR_PRE_A11: 193d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Preconditioner for the Schur complement formed from A11\n")); 194d71ae5a4SJacob Faibussowitsch break; 195d71ae5a4SJacob Faibussowitsch case PC_FIELDSPLIT_SCHUR_PRE_FULL: 196d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Preconditioner for the Schur complement formed from the exact Schur complement\n")); 197d71ae5a4SJacob Faibussowitsch break; 1983e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1993e8b8b31SMatthew G Knepley if (jac->schur_user) { 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Preconditioner for the Schur complement formed from user provided matrix\n")); 2013e8b8b31SMatthew G Knepley } else { 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Preconditioner for the Schur complement formed from A11\n")); 2033e8b8b31SMatthew G Knepley } 2043e8b8b31SMatthew G Knepley break; 205d71ae5a4SJacob Faibussowitsch default: 206d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 2073e8b8b31SMatthew G Knepley } 2089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Split info:\n")); 2099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2103b224e63SBarry Smith for (i = 0; i < jac->nsplits; i++) { 2113b224e63SBarry Smith if (ilink->fields) { 21263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Fields ", i)); 2139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 2143b224e63SBarry Smith for (j = 0; j < ilink->nfields; j++) { 21548a46eb9SPierre Jolivet if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ",")); 21663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, ilink->fields[j])); 2173b224e63SBarry Smith } 2189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 2199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 2203b224e63SBarry Smith } else { 22163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Split number %" PetscInt_FMT " Defined by IS\n", i)); 2223b224e63SBarry Smith } 2233b224e63SBarry Smith ilink = ilink->next; 2243b224e63SBarry Smith } 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for A00 block\n")); 2269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 227ac530a7eSPierre Jolivet if (jac->head) PetscCall(KSPView(jac->head->ksp, viewer)); 228ac530a7eSPierre Jolivet else PetscCall(PetscViewerASCIIPrintf(viewer, " not yet available\n")); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 23006de4afeSJed Brown if (jac->head && jac->kspupper != jac->head->ksp) { 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for upper A00 in upper triangular factor\n")); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2339566063dSJacob Faibussowitsch if (jac->kspupper) PetscCall(KSPView(jac->kspupper, viewer)); 2349566063dSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, " not yet available\n")); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 236443836d0SMatthew G Knepley } 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for S = A11 - A10 inv(A00) A01\n")); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 23912cae6f2SJed Brown if (jac->kspschur) { 2409566063dSJacob Faibussowitsch PetscCall(KSPView(jac->kspschur, viewer)); 24112cae6f2SJed Brown } else { 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " not yet available\n")); 24312cae6f2SJed Brown } 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 24606de4afeSJed Brown } else if (isdraw && jac->head) { 2474996c5bdSBarry Smith PetscDraw draw; 2484996c5bdSBarry Smith PetscReal x, y, w, wd, h; 2494996c5bdSBarry Smith PetscInt cnt = 2; 2504996c5bdSBarry Smith char str[32]; 2514996c5bdSBarry Smith 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 2539566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 254c74581afSBarry Smith if (jac->kspupper != jac->head->ksp) cnt++; 255c74581afSBarry Smith w = 2 * PetscMin(1.0 - x, x); 256c74581afSBarry Smith wd = w / (cnt + 1); 257c74581afSBarry Smith 2589566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 32, "Schur fact. %s", PCFieldSplitSchurFactTypes[jac->schurfactorization])); 2599566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 2604996c5bdSBarry Smith y -= h; 2614996c5bdSBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) { 2629566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 32, "Prec. for Schur from %s", PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11])); 2633b224e63SBarry Smith } else { 2649566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 32, "Prec. for Schur from %s", PCFieldSplitSchurPreTypes[jac->schurpre])); 2654996c5bdSBarry Smith } 2669566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x + wd * (cnt - 1) / 2.0, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 2674996c5bdSBarry Smith y -= h; 2684996c5bdSBarry Smith x = x - wd * (cnt - 1) / 2.0; 2694996c5bdSBarry Smith 2709566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, y)); 2719566063dSJacob Faibussowitsch PetscCall(KSPView(jac->head->ksp, viewer)); 2729566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 2734996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) { 2744996c5bdSBarry Smith x += wd; 2759566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, y)); 2769566063dSJacob Faibussowitsch PetscCall(KSPView(jac->kspupper, viewer)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 2784996c5bdSBarry Smith } 2794996c5bdSBarry Smith x += wd; 2809566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, y)); 2819566063dSJacob Faibussowitsch PetscCall(KSPView(jac->kspschur, viewer)); 2829566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 2833b224e63SBarry Smith } 2843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2853b224e63SBarry Smith } 2863b224e63SBarry Smith 287d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_FieldSplit_GKB(PC pc, PetscViewer viewer) 288d71ae5a4SJacob Faibussowitsch { 289a51937d4SCarola Kruse PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2909f196a02SMartin Diehl PetscBool isascii, isdraw; 291a51937d4SCarola Kruse PetscInt i, j; 292a51937d4SCarola Kruse PC_FieldSplitLink ilink = jac->head; 293a51937d4SCarola Kruse 294a51937d4SCarola Kruse PetscFunctionBegin; 2959f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 2969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 2979f196a02SMartin Diehl if (isascii) { 298a51937d4SCarola Kruse if (jac->bs > 0) { 29963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FieldSplit with %s composition: total splits = %" PetscInt_FMT ", blocksize = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits, jac->bs)); 300a51937d4SCarola Kruse } else { 30163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FieldSplit with %s composition: total splits = %" PetscInt_FMT "\n", PCCompositeTypes[jac->type], jac->nsplits)); 302a51937d4SCarola Kruse } 30348a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for blocks\n")); 30448a46eb9SPierre Jolivet if (jac->diag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for diagonal blocks\n")); 30548a46eb9SPierre Jolivet if (jac->offdiag_use_amat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat (not Pmat) as operator for off-diagonal blocks\n")); 306a51937d4SCarola Kruse 30763a3b9bcSJacob 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)); 3089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Solver info for H = A00 + nu*A01*A01' matrix:\n")); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 310a51937d4SCarola Kruse 311a51937d4SCarola Kruse if (ilink->fields) { 31263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Split number 0 Fields ")); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 314a51937d4SCarola Kruse for (j = 0; j < ilink->nfields; j++) { 31548a46eb9SPierre Jolivet if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ",")); 31663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, ilink->fields[j])); 317a51937d4SCarola Kruse } 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 320a51937d4SCarola Kruse } else { 32163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Split number 0 Defined by IS\n")); 322a51937d4SCarola Kruse } 3239566063dSJacob Faibussowitsch PetscCall(KSPView(ilink->ksp, viewer)); 324a51937d4SCarola Kruse 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 326a51937d4SCarola Kruse } 327a51937d4SCarola Kruse 328a51937d4SCarola Kruse if (isdraw) { 329a51937d4SCarola Kruse PetscDraw draw; 330a51937d4SCarola Kruse PetscReal x, y, w, wd; 331a51937d4SCarola Kruse 3329566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 3339566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 334a51937d4SCarola Kruse w = 2 * PetscMin(1.0 - x, x); 335a51937d4SCarola Kruse wd = w / (jac->nsplits + 1); 336a51937d4SCarola Kruse x = x - wd * (jac->nsplits - 1) / 2.0; 337a51937d4SCarola Kruse for (i = 0; i < jac->nsplits; i++) { 3389566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, y)); 3399566063dSJacob Faibussowitsch PetscCall(KSPView(ilink->ksp, viewer)); 3409566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 341a51937d4SCarola Kruse x += wd; 342a51937d4SCarola Kruse ilink = ilink->next; 343a51937d4SCarola Kruse } 344a51937d4SCarola Kruse } 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346a51937d4SCarola Kruse } 347a51937d4SCarola Kruse 34880670ca5SBarry Smith /* Precondition: jac->bs is set to a meaningful value or MATNEST */ 349d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 350d71ae5a4SJacob Faibussowitsch { 3516c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 35280670ca5SBarry Smith PetscInt bs, i, nfields, *ifields, nfields_col, *ifields_col; 35380670ca5SBarry Smith PetscBool flg, flg_col, mnest; 3545d4c12cdSJungho Lee char optionname[128], splitname[8], optionname_col[128]; 3556c924f48SJed Brown 3566c924f48SJed Brown PetscFunctionBegin; 35780670ca5SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat, MATNEST, &mnest)); 35880670ca5SBarry Smith if (mnest) { 35980670ca5SBarry Smith PetscCall(MatNestGetSize(pc->pmat, &bs, NULL)); 36080670ca5SBarry Smith } else { 36180670ca5SBarry Smith bs = jac->bs; 36280670ca5SBarry Smith } 36380670ca5SBarry Smith PetscCall(PetscMalloc2(bs, &ifields, bs, &ifields_col)); 3646c924f48SJed Brown for (i = 0, flg = PETSC_TRUE;; i++) { 36563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(splitname, sizeof(splitname), "%" PetscInt_FMT, i)); 36663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%" PetscInt_FMT "_fields", i)); 36763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(optionname_col, sizeof(optionname_col), "-pc_fieldsplit_%" PetscInt_FMT "_fields_col", i)); 36880670ca5SBarry Smith nfields = bs; 36980670ca5SBarry Smith nfields_col = bs; 3709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg)); 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, optionname_col, ifields_col, &nfields_col, &flg_col)); 3726c924f48SJed Brown if (!flg) break; 3735d4c12cdSJungho Lee else if (flg && !flg_col) { 37428b400f6SJacob Faibussowitsch PetscCheck(nfields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot list zero fields"); 3759566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetFields(pc, splitname, nfields, ifields, ifields)); 3762fa5cd67SKarl Rupp } else { 3777827d75bSBarry Smith PetscCheck(nfields && nfields_col, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot list zero fields"); 37808401ef6SPierre Jolivet PetscCheck(nfields == nfields_col, PETSC_COMM_SELF, PETSC_ERR_USER, "Number of row and column fields must match"); 3799566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetFields(pc, splitname, nfields, ifields, ifields_col)); 3805d4c12cdSJungho Lee } 3816c924f48SJed Brown } 3826c924f48SJed Brown if (i > 0) { 3836c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 3846c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 3856c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 3866c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 3876c924f48SJed Brown } 38880670ca5SBarry Smith PetscCall(PetscFree2(ifields, ifields_col)); 3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3906c924f48SJed Brown } 3916c924f48SJed Brown 392d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 393d71ae5a4SJacob Faibussowitsch { 3940971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 3955a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3967b752e3dSPatrick Sanan PetscBool fieldsplit_default = PETSC_FALSE, coupling = PETSC_FALSE; 3976c924f48SJed Brown PetscInt i; 3980971522cSBarry Smith 3990971522cSBarry Smith PetscFunctionBegin; 4007287d2a3SDmitry Karpeev /* 401f5f0d762SBarry Smith Kinda messy, but at least this now uses DMCreateFieldDecomposition(). 4027287d2a3SDmitry Karpeev Should probably be rewritten. 4037287d2a3SDmitry Karpeev */ 404f5f0d762SBarry Smith if (!ilink) { 4059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_fieldsplit_detect_coupling", &coupling, NULL)); 4067b752e3dSPatrick Sanan if (pc->dm && jac->dm_splits && !jac->detect && !coupling) { 407bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 4080784a22cSJed Brown char **fieldNames; 4097b62db95SJungho Lee IS *fields; 410e7c4fc90SDmitry Karpeev DM *dms; 411bafc1b83SMatthew G Knepley DM subdm[128]; 412bafc1b83SMatthew G Knepley PetscBool flg; 413bafc1b83SMatthew G Knepley 4149566063dSJacob Faibussowitsch PetscCall(DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms)); 415bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 416bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE;; i++) { 417bafc1b83SMatthew G Knepley PetscInt ifields[128]; 418bafc1b83SMatthew G Knepley IS compField; 419bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 420bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 421bafc1b83SMatthew G Knepley 42263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%" PetscInt_FMT "_fields", i)); 4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg)); 424bafc1b83SMatthew G Knepley if (!flg) break; 42563a3b9bcSJacob Faibussowitsch PetscCheck(numFields <= 128, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot currently support %" PetscInt_FMT " > 128 fields", numFields); 4269566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i])); 427bafc1b83SMatthew G Knepley if (nfields == 1) { 4289566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField)); 429bafc1b83SMatthew G Knepley } else { 43063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(splitname, sizeof(splitname), "%" PetscInt_FMT, i)); 4319566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, splitname, compField)); 4327287d2a3SDmitry Karpeev } 4339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compField)); 434bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 435bafc1b83SMatthew G Knepley f = ifields[j]; 4369566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldNames[f])); 4379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fields[f])); 4387b62db95SJungho Lee } 439bafc1b83SMatthew G Knepley } 440bafc1b83SMatthew G Knepley if (i == 0) { 441bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 4429566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, fieldNames[f], fields[f])); 4439566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldNames[f])); 4449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fields[f])); 445bafc1b83SMatthew G Knepley } 446bafc1b83SMatthew G Knepley } else { 44748a46eb9SPierre Jolivet for (j = 0; j < numFields; j++) PetscCall(DMDestroy(dms + j)); 4489566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 4499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i, &dms)); 4502fa5cd67SKarl Rupp for (j = 0; j < i; ++j) dms[j] = subdm[j]; 451bafc1b83SMatthew G Knepley } 4529566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldNames)); 4539566063dSJacob Faibussowitsch PetscCall(PetscFree(fields)); 454e7c4fc90SDmitry Karpeev if (dms) { 4559566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n")); 456bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 4577287d2a3SDmitry Karpeev const char *prefix; 458f4f49eeaSPierre Jolivet PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ilink->ksp, &prefix)); 459f4f49eeaSPierre Jolivet PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dms[i], prefix)); 4609566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ilink->ksp, dms[i])); 4619566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(ilink->ksp, PETSC_FALSE)); 4629566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)dms[i], (PetscObject)ilink->ksp, 0)); 4639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms[i])); 4642fa5ba8aSJed Brown } 4659566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 4668b8307b2SJed Brown } 46766ffff09SJed Brown } else { 468521d7252SBarry Smith if (jac->bs <= 0) { 469ac530a7eSPierre Jolivet if (pc->pmat) PetscCall(MatGetBlockSize(pc->pmat, &jac->bs)); 470ac530a7eSPierre Jolivet else jac->bs = 1; 471521d7252SBarry Smith } 472d32f9abdSBarry Smith 4737b752e3dSPatrick Sanan if (jac->detect) { 4746ce1633cSBarry Smith IS zerodiags, rest; 4756ce1633cSBarry Smith PetscInt nmin, nmax; 4766ce1633cSBarry Smith 4779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->mat, &nmin, &nmax)); 4787199da05SBarry Smith if (jac->diag_use_amat) { 4799566063dSJacob Faibussowitsch PetscCall(MatFindZeroDiagonals(pc->mat, &zerodiags)); 4807199da05SBarry Smith } else { 4819566063dSJacob Faibussowitsch PetscCall(MatFindZeroDiagonals(pc->pmat, &zerodiags)); 4827199da05SBarry Smith } 4839566063dSJacob Faibussowitsch PetscCall(ISComplement(zerodiags, nmin, nmax, &rest)); 4849566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "0", rest)); 4859566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "1", zerodiags)); 4869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&zerodiags)); 4879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rest)); 4883a062f41SBarry Smith } else if (coupling) { 4893a062f41SBarry Smith IS coupling, rest; 4903a062f41SBarry Smith PetscInt nmin, nmax; 4913a062f41SBarry Smith 4929566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->mat, &nmin, &nmax)); 4937199da05SBarry Smith if (jac->offdiag_use_amat) { 4949566063dSJacob Faibussowitsch PetscCall(MatFindOffBlockDiagonalEntries(pc->mat, &coupling)); 4957199da05SBarry Smith } else { 4969566063dSJacob Faibussowitsch PetscCall(MatFindOffBlockDiagonalEntries(pc->pmat, &coupling)); 4977199da05SBarry Smith } 4989566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc->mat), nmax - nmin, nmin, 1, &rest)); 4999566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(rest)); 5009566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "0", rest)); 5019566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "1", coupling)); 5029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coupling)); 5039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rest)); 5046ce1633cSBarry Smith } else { 5059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_fieldsplit_default", &fieldsplit_default, NULL)); 5067287d2a3SDmitry Karpeev if (!fieldsplit_default) { 507d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 508d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 5099566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetRuntimeSplits_Private(pc)); 5109566063dSJacob Faibussowitsch if (jac->splitdefined) PetscCall(PetscInfo(pc, "Splits defined using the options database\n")); 511d32f9abdSBarry Smith } 5126dbb499eSCian Wilson if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) { 5139f001fe8SStefano Zampini Mat M = pc->pmat; 514f3b928b9SStefano Zampini PetscBool isnest; 51580670ca5SBarry Smith PetscInt nf; 516f3b928b9SStefano Zampini 5179566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Using default splitting of fields\n")); 5189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATNEST, &isnest)); 519f3b928b9SStefano Zampini if (!isnest) { 5209f001fe8SStefano Zampini M = pc->mat; 5219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat, MATNEST, &isnest)); 522f3b928b9SStefano Zampini } 52380670ca5SBarry Smith if (!isnest) nf = jac->bs; 52480670ca5SBarry Smith else PetscCall(MatNestGetSize(M, &nf, NULL)); 52580670ca5SBarry Smith for (i = 0; i < nf; i++) { 5266c924f48SJed Brown char splitname[8]; 52780670ca5SBarry Smith 52863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(splitname, sizeof(splitname), "%" PetscInt_FMT, i)); 5299566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetFields(pc, splitname, 1, &i, &i)); 53079416396SBarry Smith } 5315d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 532521d7252SBarry Smith } 53366ffff09SJed Brown } 5346ce1633cSBarry Smith } 535edf189efSBarry Smith } else if (jac->nsplits == 1) { 536edf189efSBarry Smith IS is2; 537edf189efSBarry Smith PetscInt nmin, nmax; 538edf189efSBarry Smith 5390fdf79fbSJacob Faibussowitsch PetscCheck(ilink->is, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Must provide at least two sets of fields to PCFieldSplit()"); 5409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->mat, &nmin, &nmax)); 5419566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is, nmin, nmax, &is2)); 5429566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "1", is2)); 5439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 544edf189efSBarry Smith } 545d0af7cd3SBarry Smith 54663a3b9bcSJacob Faibussowitsch PetscCheck(jac->nsplits >= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Unhandled case, must have at least two fields, not %" PetscInt_FMT, jac->nsplits); 5473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54869a612a9SBarry Smith } 54969a612a9SBarry Smith 550d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGolubKahanComputeExplicitOperator(Mat A, Mat B, Mat C, Mat *H, PetscReal gkbnu) 551d71ae5a4SJacob Faibussowitsch { 552a51937d4SCarola Kruse Mat BT, T; 553de482cd7SCarola Kruse PetscReal nrmT, nrmB; 554a51937d4SCarola Kruse 555a51937d4SCarola Kruse PetscFunctionBegin; 5569566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(C, MAT_INITIAL_MATRIX, &T)); /* Test if augmented matrix is symmetric */ 5579566063dSJacob Faibussowitsch PetscCall(MatAXPY(T, -1.0, B, DIFFERENT_NONZERO_PATTERN)); 5589566063dSJacob Faibussowitsch PetscCall(MatNorm(T, NORM_1, &nrmT)); 5599566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_1, &nrmB)); 560*b0c98d1dSPierre Jolivet PetscCheck(nrmB <= 0 || nrmT / nrmB < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix is not symmetric/Hermitian, GKB is not applicable."); 561049d1499SBarry Smith 562a51937d4SCarola Kruse /* Compute augmented Lagrangian matrix H = A00 + nu*A01*A01'. This corresponds to */ 563a51937d4SCarola Kruse /* setting N := 1/nu*I in [Ar13]. */ 5649566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_INITIAL_MATRIX, &BT)); 565fb842aefSJose E. Roman PetscCall(MatMatMult(B, BT, MAT_INITIAL_MATRIX, PETSC_CURRENT, H)); /* H = A01*A01' */ 5669566063dSJacob Faibussowitsch PetscCall(MatAYPX(*H, gkbnu, A, DIFFERENT_NONZERO_PATTERN)); /* H = A00 + nu*A01*A01' */ 567a51937d4SCarola Kruse 5689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 5699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 571a51937d4SCarola Kruse } 572a51937d4SCarola Kruse 57354a546c1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions, const char pre[], const char name[], const char *option[], const char *value[], PetscBool *flg); 574514bf10dSMatthew G Knepley 575d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_FieldSplit(PC pc) 576d71ae5a4SJacob Faibussowitsch { 57769a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 5785a9f2f41SSatish Balay PC_FieldSplitLink ilink; 5792c9966d7SBarry Smith PetscInt i, nsplit; 5802033cbf1SStefano Zampini PetscBool matnest = PETSC_FALSE; 58169a612a9SBarry Smith 58269a612a9SBarry Smith PetscFunctionBegin; 5835da88fe4STristan Konolige pc->failedreason = PC_NOERROR; 5849566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetDefaults(pc)); 58597bbdb24SBarry Smith nsplit = jac->nsplits; 5865a9f2f41SSatish Balay ilink = jac->head; 58780670ca5SBarry Smith if (pc->pmat) PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATNEST, &matnest)); 58897bbdb24SBarry Smith 58997bbdb24SBarry Smith /* get the matrices for each split */ 590704ba839SBarry Smith if (!jac->issetup) { 5911b9fc7fcSBarry Smith PetscInt rstart, rend, nslots, bs; 59297bbdb24SBarry Smith 593704ba839SBarry Smith jac->issetup = PETSC_TRUE; 594704ba839SBarry Smith 5955d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 5962c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 5972c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 5982c9966d7SBarry Smith } 5994db63379SBarry Smith 6004db63379SBarry Smith /* MatCreateSubMatrix() for [S]BAIJ matrices can only work if the indices include entire blocks of the matrix */ 6014db63379SBarry Smith PetscCall(MatGetBlockSize(pc->pmat, &bs)); 6024db63379SBarry Smith if (bs > 1 && (jac->bs <= bs || jac->bs % bs)) { 6034db63379SBarry Smith PetscBool blk; 6044db63379SBarry Smith 6054db63379SBarry Smith PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &blk, MATBAIJ, MATSBAIJ, MATSEQBAIJ, MATSEQSBAIJ, MATMPIBAIJ, MATMPISBAIJ, NULL)); 6064db63379SBarry Smith PetscCheck(!blk, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Cannot use MATBAIJ with PCFIELDSPLIT and currently set matrix and PC blocksizes"); 6074db63379SBarry Smith } 6084db63379SBarry Smith 60980670ca5SBarry Smith if (!matnest) { /* use the matrix blocksize and stride IS to determine the index sets that define the submatrices */ 61051f519a2SBarry Smith bs = jac->bs; 6119566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat, &rstart, &rend)); 6121b9fc7fcSBarry Smith nslots = (rend - rstart) / bs; 6131b9fc7fcSBarry Smith for (i = 0; i < nsplit; i++) { 6141b9fc7fcSBarry Smith if (jac->defaultsplit) { 6159566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), nslots, rstart + i, nsplit, &ilink->is)); 6162033cbf1SStefano Zampini PetscCall(PetscObjectReference((PetscObject)ilink->is)); 6172033cbf1SStefano Zampini ilink->is_col = ilink->is; 618704ba839SBarry Smith } else if (!ilink->is) { 6192033cbf1SStefano Zampini PetscBool same_fields = PETSC_TRUE; 6202033cbf1SStefano Zampini 6212033cbf1SStefano Zampini for (PetscInt k = 0; k < ilink->nfields; k++) { 6222033cbf1SStefano Zampini if (ilink->fields[k] != ilink->fields_col[k]) same_fields = PETSC_FALSE; 6232033cbf1SStefano Zampini } 6242033cbf1SStefano Zampini 625ccb205f8SBarry Smith if (ilink->nfields > 1) { 6265f4ab4e1SJungho Lee PetscInt *ii, *jj, j, k, nfields = ilink->nfields, *fields = ilink->fields, *fields_col = ilink->fields_col; 62780670ca5SBarry Smith 6289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ilink->nfields * nslots, &ii)); 6292033cbf1SStefano Zampini if (!same_fields) PetscCall(PetscMalloc1(ilink->nfields * nslots, &jj)); 6301b9fc7fcSBarry Smith for (j = 0; j < nslots; j++) { 6311b9fc7fcSBarry Smith for (k = 0; k < nfields; k++) { 6321b9fc7fcSBarry Smith ii[nfields * j + k] = rstart + bs * j + fields[k]; 6332033cbf1SStefano Zampini if (!same_fields) jj[nfields * j + k] = rstart + bs * j + fields_col[k]; 63497bbdb24SBarry Smith } 63597bbdb24SBarry Smith } 6369566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), nslots * nfields, ii, PETSC_OWN_POINTER, &ilink->is)); 6372033cbf1SStefano Zampini if (!same_fields) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), nslots * nfields, jj, PETSC_OWN_POINTER, &ilink->is_col)); 6382033cbf1SStefano Zampini else { 6392033cbf1SStefano Zampini PetscCall(PetscObjectReference((PetscObject)ilink->is)); 6402033cbf1SStefano Zampini ilink->is_col = ilink->is; 6412033cbf1SStefano Zampini } 6429566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(ilink->is, nfields)); 6439566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(ilink->is_col, nfields)); 644ccb205f8SBarry Smith } else { 6459566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), nslots, rstart + ilink->fields[0], bs, &ilink->is)); 6462033cbf1SStefano Zampini if (!same_fields) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), nslots, rstart + ilink->fields_col[0], bs, &ilink->is_col)); 6472033cbf1SStefano Zampini else { 6482033cbf1SStefano Zampini PetscCall(PetscObjectReference((PetscObject)ilink->is)); 6492033cbf1SStefano Zampini ilink->is_col = ilink->is; 650ccb205f8SBarry Smith } 6513e197d65SBarry Smith } 6522033cbf1SStefano Zampini } 653704ba839SBarry Smith ilink = ilink->next; 6541b9fc7fcSBarry Smith } 65580670ca5SBarry Smith } else { /* use the IS that define the MATNEST to determine the index sets that define the submatrices */ 65680670ca5SBarry Smith IS *rowis, *colis, *ises = NULL; 65780670ca5SBarry Smith PetscInt mis, nis; 65880670ca5SBarry Smith 65980670ca5SBarry Smith PetscCall(MatNestGetSize(pc->pmat, &mis, &nis)); 66080670ca5SBarry Smith PetscCall(PetscMalloc2(mis, &rowis, nis, &colis)); 66180670ca5SBarry Smith PetscCall(MatNestGetISs(pc->pmat, rowis, colis)); 66280670ca5SBarry Smith if (!jac->defaultsplit) PetscCall(PetscMalloc1(mis, &ises)); 66380670ca5SBarry Smith 66480670ca5SBarry Smith for (i = 0; i < nsplit; i++) { 66580670ca5SBarry Smith if (jac->defaultsplit) { 66680670ca5SBarry Smith PetscCall(ISDuplicate(rowis[i], &ilink->is)); 6672033cbf1SStefano Zampini PetscCall(PetscObjectReference((PetscObject)ilink->is)); 6682033cbf1SStefano Zampini ilink->is_col = ilink->is; 66980670ca5SBarry Smith } else if (!ilink->is) { 67080670ca5SBarry Smith if (ilink->nfields > 1) { 67180670ca5SBarry Smith for (PetscInt j = 0; j < ilink->nfields; j++) ises[j] = rowis[ilink->fields[j]]; 67280670ca5SBarry Smith PetscCall(ISConcatenate(PetscObjectComm((PetscObject)pc), ilink->nfields, ises, &ilink->is)); 67380670ca5SBarry Smith } else { 67480670ca5SBarry Smith PetscCall(ISDuplicate(rowis[ilink->fields[0]], &ilink->is)); 67580670ca5SBarry Smith } 6762033cbf1SStefano Zampini PetscCall(PetscObjectReference((PetscObject)ilink->is)); 6772033cbf1SStefano Zampini ilink->is_col = ilink->is; 67880670ca5SBarry Smith } 67980670ca5SBarry Smith ilink = ilink->next; 68080670ca5SBarry Smith } 68180670ca5SBarry Smith PetscCall(PetscFree2(rowis, colis)); 68280670ca5SBarry Smith PetscCall(PetscFree(ises)); 68380670ca5SBarry Smith } 6841b9fc7fcSBarry Smith } 6851b9fc7fcSBarry Smith 686704ba839SBarry Smith ilink = jac->head; 68797bbdb24SBarry Smith if (!jac->pmat) { 688bdddcaaaSMatthew G Knepley Vec xtmp; 689bdddcaaaSMatthew G Knepley 6909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &xtmp, NULL)); 6919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsplit, &jac->pmat)); 6929566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nsplit, &jac->x, nsplit, &jac->y)); 693cf502942SBarry Smith for (i = 0; i < nsplit; i++) { 694bdddcaaaSMatthew G Knepley MatNullSpace sp; 695bdddcaaaSMatthew G Knepley 6967addb90fSBarry Smith /* Check for matrix attached to IS */ 6979566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ilink->is, "pmat", (PetscObject *)&jac->pmat[i])); 698a3df900dSMatthew G Knepley if (jac->pmat[i]) { 6999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)jac->pmat[i])); 700a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 701a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 7022fa5cd67SKarl Rupp 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)jac->schur_user)); 704a3df900dSMatthew G Knepley } 705a3df900dSMatthew G Knepley } else { 7063a062f41SBarry Smith const char *prefix; 7079566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ilink->is_col, MAT_INITIAL_MATRIX, &jac->pmat[i])); 7082f427464SPierre Jolivet PetscCall(MatGetOptionsPrefix(jac->pmat[i], &prefix)); 7092f427464SPierre Jolivet if (!prefix) { 7109566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(ilink->ksp, &prefix)); 7119566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(jac->pmat[i], prefix)); 7122f427464SPierre Jolivet } 71345881c45SPierre Jolivet PetscCall(MatSetFromOptions(jac->pmat[i])); 7149566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(jac->pmat[i], NULL, "-mat_view")); 715a3df900dSMatthew G Knepley } 716bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 7179566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(jac->pmat[i], &jac->x[i], &jac->y[i])); 7189371c9d4SSatish Balay ilink->x = jac->x[i]; 7199371c9d4SSatish Balay ilink->y = jac->y[i]; 7209371c9d4SSatish Balay ilink->z = NULL; 721bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 7229566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(xtmp, ilink->is, jac->x[i], NULL, &ilink->sctx)); 7239566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ilink->is, "nearnullspace", (PetscObject *)&sp)); 72448a46eb9SPierre Jolivet if (sp) PetscCall(MatSetNearNullSpace(jac->pmat[i], sp)); 725704ba839SBarry Smith ilink = ilink->next; 726cf502942SBarry Smith } 7279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xtmp)); 72897bbdb24SBarry Smith } else { 729ef7efd37SHong Zhang MatReuse scall; 7304849c82aSBarry Smith MatNullSpace *nullsp = NULL; 7314849c82aSBarry Smith 732ef7efd37SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 7334849c82aSBarry Smith PetscCall(MatGetNullSpaces(nsplit, jac->pmat, &nullsp)); 73448a46eb9SPierre Jolivet for (i = 0; i < nsplit; i++) PetscCall(MatDestroy(&jac->pmat[i])); 735ef7efd37SHong Zhang scall = MAT_INITIAL_MATRIX; 736ef7efd37SHong Zhang } else scall = MAT_REUSE_MATRIX; 737ef7efd37SHong Zhang 738cf502942SBarry Smith for (i = 0; i < nsplit; i++) { 739a3df900dSMatthew G Knepley Mat pmat; 740a3df900dSMatthew G Knepley 7417addb90fSBarry Smith /* Check for matrix attached to IS */ 7429566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ilink->is, "pmat", (PetscObject *)&pmat)); 74348a46eb9SPierre Jolivet if (!pmat) PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ilink->is_col, scall, &jac->pmat[i])); 744704ba839SBarry Smith ilink = ilink->next; 745cf502942SBarry Smith } 7464849c82aSBarry Smith if (nullsp) PetscCall(MatRestoreNullSpaces(nsplit, jac->pmat, &nullsp)); 74797bbdb24SBarry Smith } 7484e39094bSDmitry Karpeev if (jac->diag_use_amat) { 749519d70e2SJed Brown ilink = jac->head; 750519d70e2SJed Brown if (!jac->mat) { 7519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsplit, &jac->mat)); 752519d70e2SJed Brown for (i = 0; i < nsplit; i++) { 7539566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ilink->is_col, MAT_INITIAL_MATRIX, &jac->mat[i])); 754519d70e2SJed Brown ilink = ilink->next; 755519d70e2SJed Brown } 756519d70e2SJed Brown } else { 757ef7efd37SHong Zhang MatReuse scall; 7584849c82aSBarry Smith MatNullSpace *nullsp = NULL; 7594849c82aSBarry Smith 760ef7efd37SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 7614849c82aSBarry Smith PetscCall(MatGetNullSpaces(nsplit, jac->mat, &nullsp)); 76248a46eb9SPierre Jolivet for (i = 0; i < nsplit; i++) PetscCall(MatDestroy(&jac->mat[i])); 763ef7efd37SHong Zhang scall = MAT_INITIAL_MATRIX; 764ef7efd37SHong Zhang } else scall = MAT_REUSE_MATRIX; 765ef7efd37SHong Zhang 766ef7efd37SHong Zhang for (i = 0; i < nsplit; i++) { 7679566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ilink->is_col, scall, &jac->mat[i])); 768519d70e2SJed Brown ilink = ilink->next; 769519d70e2SJed Brown } 7704849c82aSBarry Smith if (nullsp) PetscCall(MatRestoreNullSpaces(nsplit, jac->mat, &nullsp)); 771519d70e2SJed Brown } 772519d70e2SJed Brown } else { 773519d70e2SJed Brown jac->mat = jac->pmat; 774519d70e2SJed Brown } 77597bbdb24SBarry Smith 77653935eafSBarry Smith /* Check for null space attached to IS */ 77753935eafSBarry Smith ilink = jac->head; 77853935eafSBarry Smith for (i = 0; i < nsplit; i++) { 77953935eafSBarry Smith MatNullSpace sp; 78053935eafSBarry Smith 7819566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ilink->is, "nullspace", (PetscObject *)&sp)); 78248a46eb9SPierre Jolivet if (sp) PetscCall(MatSetNullSpace(jac->mat[i], sp)); 78353935eafSBarry Smith ilink = ilink->next; 78453935eafSBarry Smith } 78553935eafSBarry Smith 786a51937d4SCarola Kruse if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR && jac->type != PC_COMPOSITE_GKB) { 78768dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 7884e39094bSDmitry Karpeev /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */ 78968dd23aaSBarry Smith ilink = jac->head; 790e52d2c62SBarry Smith if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) { 791e52d2c62SBarry 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 */ 792e52d2c62SBarry Smith if (!jac->Afield) { 7939566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nsplit, &jac->Afield)); 79480c96bb1SFande Kong if (jac->offdiag_use_amat) { 7959566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->next->is, ilink->is, MAT_INITIAL_MATRIX, &jac->Afield[1])); 796e52d2c62SBarry Smith } else { 7979566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->next->is, ilink->is, MAT_INITIAL_MATRIX, &jac->Afield[1])); 79880c96bb1SFande Kong } 79980c96bb1SFande Kong } else { 800ef7efd37SHong Zhang MatReuse scall; 801e9422dd5SStefano Zampini 802ef7efd37SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 8039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->Afield[1])); 804ef7efd37SHong Zhang scall = MAT_INITIAL_MATRIX; 805ef7efd37SHong Zhang } else scall = MAT_REUSE_MATRIX; 806ef7efd37SHong Zhang 80780c96bb1SFande Kong if (jac->offdiag_use_amat) { 8089566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->next->is, ilink->is, scall, &jac->Afield[1])); 80980c96bb1SFande Kong } else { 8109566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->next->is, ilink->is, scall, &jac->Afield[1])); 81180c96bb1SFande Kong } 812e52d2c62SBarry Smith } 813e52d2c62SBarry Smith } else { 81468dd23aaSBarry Smith if (!jac->Afield) { 8159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsplit, &jac->Afield)); 81668dd23aaSBarry Smith for (i = 0; i < nsplit; i++) { 81780c96bb1SFande Kong if (jac->offdiag_use_amat) { 8189566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, NULL, MAT_INITIAL_MATRIX, &jac->Afield[i])); 81980c96bb1SFande Kong } else { 8209566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, NULL, MAT_INITIAL_MATRIX, &jac->Afield[i])); 82180c96bb1SFande Kong } 82268dd23aaSBarry Smith ilink = ilink->next; 82368dd23aaSBarry Smith } 82468dd23aaSBarry Smith } else { 825ef7efd37SHong Zhang MatReuse scall; 826ef7efd37SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 82748a46eb9SPierre Jolivet for (i = 0; i < nsplit; i++) PetscCall(MatDestroy(&jac->Afield[i])); 828ef7efd37SHong Zhang scall = MAT_INITIAL_MATRIX; 829ef7efd37SHong Zhang } else scall = MAT_REUSE_MATRIX; 830ef7efd37SHong Zhang 83168dd23aaSBarry Smith for (i = 0; i < nsplit; i++) { 83280c96bb1SFande Kong if (jac->offdiag_use_amat) { 8339566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, NULL, scall, &jac->Afield[i])); 83480c96bb1SFande Kong } else { 8359566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, NULL, scall, &jac->Afield[i])); 83680c96bb1SFande Kong } 83768dd23aaSBarry Smith ilink = ilink->next; 83868dd23aaSBarry Smith } 83968dd23aaSBarry Smith } 84068dd23aaSBarry Smith } 841e52d2c62SBarry Smith } 84268dd23aaSBarry Smith 8433b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 8443b224e63SBarry Smith IS ccis; 84518f54938SStefano Zampini PetscBool isset, isspd = PETSC_FALSE, issym = PETSC_FALSE, flg; 8464aa3045dSJed Brown PetscInt rstart, rend; 847093c86ffSJed Brown char lscname[256]; 848093c86ffSJed Brown PetscObject LSC_L; 849ce94432eSBarry Smith 85008401ef6SPierre Jolivet PetscCheck(nsplit == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "To use Schur complement preconditioner you must have exactly 2 fields"); 85168dd23aaSBarry Smith 852c096484dSStefano Zampini /* If pc->mat is SPD, don't scale by -1 the Schur complement */ 853b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd)); 85418f54938SStefano Zampini if (jac->schurscale == (PetscScalar)-1.0) jac->schurscale = (isset && isspd) ? 1.0 : -1.0; 85518f54938SStefano Zampini PetscCall(MatIsSymmetricKnown(pc->pmat, &isset, &issym)); 856c096484dSStefano Zampini 857e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 8589566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(pc->mat, &rstart, &rend)); 859558f3fe8SPierre Jolivet PetscCall(PetscObjectTypeCompareAny(jac->offdiag_use_amat ? (PetscObject)pc->mat : (PetscObject)pc->pmat, &flg, MATSEQSBAIJ, MATMPISBAIJ, "")); 860e6cab6aaSJed Brown 8613b224e63SBarry Smith if (jac->schur) { 8620298fd71SBarry Smith KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper; 863e9422dd5SStefano Zampini MatReuse scall; 864e9422dd5SStefano Zampini 865e9422dd5SStefano Zampini if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 866e9422dd5SStefano Zampini scall = MAT_INITIAL_MATRIX; 8679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->B)); 8689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->C)); 869e9422dd5SStefano Zampini } else scall = MAT_REUSE_MATRIX; 870443836d0SMatthew G Knepley 8719566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(jac->schur, &kspInner)); 8723b224e63SBarry Smith ilink = jac->head; 8739566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); 8744e39094bSDmitry Karpeev if (jac->offdiag_use_amat) { 8759566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, scall, &jac->B)); 876c84da90fSDmitry Karpeev } else { 8779566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, scall, &jac->B)); 878c84da90fSDmitry Karpeev } 8799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ccis)); 880558f3fe8SPierre Jolivet if (!flg) { 8813b224e63SBarry Smith ilink = ilink->next; 8829566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); 8834e39094bSDmitry Karpeev if (jac->offdiag_use_amat) { 8849566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, scall, &jac->C)); 885c84da90fSDmitry Karpeev } else { 8869566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, scall, &jac->C)); 887c84da90fSDmitry Karpeev } 8889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ccis)); 889558f3fe8SPierre Jolivet } else { 89018f54938SStefano Zampini PetscCall(MatIsHermitianKnown(jac->offdiag_use_amat ? pc->mat : pc->pmat, &isset, &flg)); 89118f54938SStefano Zampini if (isset && flg) PetscCall(MatCreateHermitianTranspose(jac->B, &jac->C)); 892558f3fe8SPierre Jolivet else PetscCall(MatCreateTranspose(jac->B, &jac->C)); 893558f3fe8SPierre Jolivet } 8949566063dSJacob Faibussowitsch PetscCall(MatSchurComplementUpdateSubMatrices(jac->schur, jac->mat[0], jac->pmat[0], jac->B, jac->C, jac->mat[1])); 895a7476a74SDmitry Karpeev if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) { 8969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->schurp)); 8979566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetPmat(jac->schur, MAT_INITIAL_MATRIX, &jac->schurp)); 8985becce15SPierre Jolivet } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL && jac->kspupper != jac->head->ksp) { 899d9eadc85SPierre Jolivet PetscCall(MatDestroy(&jac->schur_user)); 900d9eadc85SPierre Jolivet PetscCall(MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user)); 901a7476a74SDmitry Karpeev } 90248a46eb9SPierre Jolivet if (kspA != kspInner) PetscCall(KSPSetOperators(kspA, jac->mat[0], jac->pmat[0])); 90348a46eb9SPierre Jolivet if (kspUpper != kspA) PetscCall(KSPSetOperators(kspUpper, jac->mat[0], jac->pmat[0])); 9049566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->kspschur, jac->schur, FieldSplitSchurPre(jac))); 9053b224e63SBarry Smith } else { 906bafc1b83SMatthew G Knepley const char *Dprefix; 907470b340bSDmitry Karpeev char schurprefix[256], schurmatprefix[256]; 908514bf10dSMatthew G Knepley char schurtestoption[256]; 909bdddcaaaSMatthew G Knepley MatNullSpace sp; 910686bed4dSStefano Zampini KSP kspt; 9113b224e63SBarry Smith 912a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 9133b224e63SBarry Smith ilink = jac->head; 9149566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); 9154e39094bSDmitry Karpeev if (jac->offdiag_use_amat) { 9169566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B)); 917c84da90fSDmitry Karpeev } else { 9189566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B)); 919c84da90fSDmitry Karpeev } 9209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ccis)); 9213b224e63SBarry Smith ilink = ilink->next; 922558f3fe8SPierre Jolivet if (!flg) { 9239566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); 9244e39094bSDmitry Karpeev if (jac->offdiag_use_amat) { 9259566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C)); 926c84da90fSDmitry Karpeev } else { 9279566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C)); 928c84da90fSDmitry Karpeev } 9299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ccis)); 930558f3fe8SPierre Jolivet } else { 93118f54938SStefano Zampini PetscCall(MatIsHermitianKnown(jac->offdiag_use_amat ? pc->mat : pc->pmat, &isset, &flg)); 93218f54938SStefano Zampini if (isset && flg) PetscCall(MatCreateHermitianTranspose(jac->B, &jac->C)); 933558f3fe8SPierre Jolivet else PetscCall(MatCreateTranspose(jac->B, &jac->C)); 934558f3fe8SPierre Jolivet } 935f5236f50SJed Brown /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */ 9369566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)jac->mat[0])->comm, &jac->schur)); 9379566063dSJacob Faibussowitsch PetscCall(MatSetType(jac->schur, MATSCHURCOMPLEMENT)); 9389566063dSJacob Faibussowitsch PetscCall(MatSchurComplementSetSubMatrices(jac->schur, jac->mat[0], jac->pmat[0], jac->B, jac->C, jac->mat[1])); 9399566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname)); 9409566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(jac->schur, schurmatprefix)); 9419566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(jac->schur, &kspt)); 9429566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(kspt, schurmatprefix)); 943686bed4dSStefano Zampini 944686bed4dSStefano Zampini /* Note: this is not true in general */ 9459566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(jac->mat[1], &sp)); 9461baa6e33SBarry Smith if (sp) PetscCall(MatSetNullSpace(jac->schur, sp)); 94720252d06SBarry Smith 9489566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname)); 94954a546c1SMatthew G. Knepley PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, schurtestoption, NULL, NULL, &flg)); 950514bf10dSMatthew G Knepley if (flg) { 951514bf10dSMatthew G Knepley DM dmInner; 95221635b76SJed Brown KSP kspInner; 953686bed4dSStefano Zampini PC pcInner; 95421635b76SJed Brown 9559566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(jac->schur, &kspInner)); 9569566063dSJacob Faibussowitsch PetscCall(KSPReset(kspInner)); 9579566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(kspInner, jac->mat[0], jac->pmat[0])); 9589566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname)); 95921635b76SJed Brown /* Indent this deeper to emphasize the "inner" nature of this solver. */ 9609566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject)pc, 2)); 9619566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspInner->pc, (PetscObject)pc, 2)); 9629566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(kspInner, schurprefix)); 963514bf10dSMatthew G Knepley 964514bf10dSMatthew G Knepley /* Set DM for new solver */ 9659566063dSJacob Faibussowitsch PetscCall(KSPGetDM(jac->head->ksp, &dmInner)); 9669566063dSJacob Faibussowitsch PetscCall(KSPSetDM(kspInner, dmInner)); 9679566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(kspInner, PETSC_FALSE)); 968686bed4dSStefano Zampini 969686bed4dSStefano Zampini /* Defaults to PCKSP as preconditioner */ 9709566063dSJacob Faibussowitsch PetscCall(KSPGetPC(kspInner, &pcInner)); 9719566063dSJacob Faibussowitsch PetscCall(PCSetType(pcInner, PCKSP)); 9729566063dSJacob Faibussowitsch PetscCall(PCKSPSetKSP(pcInner, jac->head->ksp)); 973514bf10dSMatthew G Knepley } else { 97421635b76SJed Brown /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or 97521635b76SJed Brown * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact, 97621635b76SJed Brown * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for 97721635b76SJed 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 97821635b76SJed Brown * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used 97921635b76SJed Brown * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */ 9809566063dSJacob Faibussowitsch PetscCall(KSPSetType(jac->head->ksp, KSPGMRES)); 9819566063dSJacob Faibussowitsch PetscCall(MatSchurComplementSetKSP(jac->schur, jac->head->ksp)); 982bafc1b83SMatthew G Knepley } 9839566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->head->ksp, jac->mat[0], jac->pmat[0])); 9849566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(jac->head->ksp)); 9859566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(jac->schur)); 9863b224e63SBarry Smith 9879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)jac->schur, MATSCHURCOMPLEMENT, &flg)); 988686bed4dSStefano Zampini if (flg) { /* Need to do this otherwise PCSetUp_KSP will overwrite the amat of jac->head->ksp */ 989686bed4dSStefano Zampini KSP kspInner; 990686bed4dSStefano Zampini PC pcInner; 991686bed4dSStefano Zampini 9929566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(jac->schur, &kspInner)); 9939566063dSJacob Faibussowitsch PetscCall(KSPGetPC(kspInner, &pcInner)); 9949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pcInner, PCKSP, &flg)); 995686bed4dSStefano Zampini if (flg) { 996686bed4dSStefano Zampini KSP ksp; 997686bed4dSStefano Zampini 9989566063dSJacob Faibussowitsch PetscCall(PCKSPGetKSP(pcInner, &ksp)); 99948a46eb9SPierre Jolivet if (ksp == jac->head->ksp) PetscCall(PCSetUseAmat(pcInner, PETSC_TRUE)); 1000686bed4dSStefano Zampini } 1001686bed4dSStefano Zampini } 10029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname)); 100354a546c1SMatthew G. Knepley PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, schurtestoption, NULL, NULL, &flg)); 1004443836d0SMatthew G Knepley if (flg) { 1005443836d0SMatthew G Knepley DM dmInner; 1006443836d0SMatthew G Knepley 10079566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname)); 10089566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper)); 10093821be0aSBarry Smith PetscCall(KSPSetNestLevel(jac->kspupper, pc->kspnestlevel)); 10109566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(jac->kspupper, pc->erroriffailure)); 10119566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(jac->kspupper, schurprefix)); 10129566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspupper, (PetscObject)pc, 1)); 10139566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspupper->pc, (PetscObject)pc, 1)); 10149566063dSJacob Faibussowitsch PetscCall(KSPGetDM(jac->head->ksp, &dmInner)); 10159566063dSJacob Faibussowitsch PetscCall(KSPSetDM(jac->kspupper, dmInner)); 10169566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(jac->kspupper, PETSC_FALSE)); 10179566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(jac->kspupper)); 10189566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->kspupper, jac->mat[0], jac->pmat[0])); 10199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(jac->head->x, &jac->head->z)); 1020443836d0SMatthew G Knepley } else { 1021443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 10229566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)jac->head->ksp)); 1023443836d0SMatthew G Knepley } 1024443836d0SMatthew G Knepley 102548a46eb9SPierre Jolivet if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) PetscCall(MatSchurComplementGetPmat(jac->schur, MAT_INITIAL_MATRIX, &jac->schurp)); 10269566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspschur)); 10273821be0aSBarry Smith PetscCall(KSPSetNestLevel(jac->kspschur, pc->kspnestlevel)); 10289566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(jac->kspschur, pc->erroriffailure)); 10299566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspschur, (PetscObject)pc, 1)); 1030084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 10317233a360SDmitry Karpeev PC pcschur; 10329566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->kspschur, &pcschur)); 10339566063dSJacob Faibussowitsch PetscCall(PCSetType(pcschur, PCNONE)); 1034084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 1035e74569cdSMatthew G. Knepley } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) { 10365becce15SPierre Jolivet if (jac->schurfactorization != PC_FIELDSPLIT_SCHUR_FACT_FULL || jac->kspupper != jac->head->ksp) PetscCall(MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user)); 1037e69d4d44SBarry Smith } 10389566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->kspschur, jac->schur, FieldSplitSchurPre(jac))); 10399566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix)); 10409566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(jac->kspschur, Dprefix)); 1041c096484dSStefano Zampini /* propagate DM */ 1042b20b4189SMatthew G. Knepley { 1043b20b4189SMatthew G. Knepley DM sdm; 10449566063dSJacob Faibussowitsch PetscCall(KSPGetDM(jac->head->next->ksp, &sdm)); 1045b20b4189SMatthew G. Knepley if (sdm) { 10469566063dSJacob Faibussowitsch PetscCall(KSPSetDM(jac->kspschur, sdm)); 10479566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(jac->kspschur, PETSC_FALSE)); 1048b20b4189SMatthew G. Knepley } 1049b20b4189SMatthew G. Knepley } 10503b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 105120b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 10529566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(jac->kspschur)); 10533b224e63SBarry Smith } 10549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac->schur, MAT_FINAL_ASSEMBLY)); 10559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac->schur, MAT_FINAL_ASSEMBLY)); 105618f54938SStefano Zampini if (issym) PetscCall(MatSetOption(jac->schur, MAT_SYMMETRIC, PETSC_TRUE)); 105718f54938SStefano Zampini if (isspd) PetscCall(MatSetOption(jac->schur, MAT_SPD, PETSC_TRUE)); 1058093c86ffSJed Brown 1059093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 10609566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(lscname, sizeof(lscname), "%s_LSC_L", ilink->splitname)); 1061835f2295SStefano Zampini PetscCall(PetscObjectQuery((PetscObject)pc->mat, lscname, &LSC_L)); 1062835f2295SStefano Zampini if (!LSC_L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, lscname, &LSC_L)); 1063835f2295SStefano Zampini if (LSC_L) PetscCall(PetscObjectCompose((PetscObject)jac->schur, "LSC_L", LSC_L)); 10649566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(lscname, sizeof(lscname), "%s_LSC_Lp", ilink->splitname)); 1065835f2295SStefano Zampini PetscCall(PetscObjectQuery((PetscObject)pc->pmat, lscname, &LSC_L)); 1066835f2295SStefano Zampini if (!LSC_L) PetscCall(PetscObjectQuery((PetscObject)pc->mat, lscname, &LSC_L)); 1067835f2295SStefano Zampini if (LSC_L) PetscCall(PetscObjectCompose((PetscObject)jac->schur, "LSC_Lp", LSC_L)); 1068a51937d4SCarola Kruse } else if (jac->type == PC_COMPOSITE_GKB) { 1069a51937d4SCarola Kruse IS ccis; 1070a51937d4SCarola Kruse PetscInt rstart, rend; 1071a51937d4SCarola Kruse 107208401ef6SPierre Jolivet PetscCheck(nsplit == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "To use GKB preconditioner you must have exactly 2 fields"); 1073a51937d4SCarola Kruse 1074a51937d4SCarola Kruse ilink = jac->head; 1075a51937d4SCarola Kruse 1076a51937d4SCarola Kruse /* When extracting off-diagonal submatrices, we take complements from this range */ 10779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(pc->mat, &rstart, &rend)); 1078a51937d4SCarola Kruse 10799566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); 1080a51937d4SCarola Kruse if (jac->offdiag_use_amat) { 10819566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B)); 1082a51937d4SCarola Kruse } else { 10839566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B)); 1084a51937d4SCarola Kruse } 10859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ccis)); 1086e071a0a4SCarola Kruse /* Create work vectors for GKB algorithm */ 10879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ilink->x, &jac->u)); 10889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ilink->x, &jac->Hu)); 10899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ilink->x, &jac->w2)); 1090a51937d4SCarola Kruse ilink = ilink->next; 10919566063dSJacob Faibussowitsch PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); 1092a51937d4SCarola Kruse if (jac->offdiag_use_amat) { 10939566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->mat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C)); 1094a51937d4SCarola Kruse } else { 10959566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->C)); 1096a51937d4SCarola Kruse } 10979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ccis)); 1098e071a0a4SCarola Kruse /* Create work vectors for GKB algorithm */ 10999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ilink->x, &jac->v)); 11009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ilink->x, &jac->d)); 11019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ilink->x, &jac->w1)); 11029566063dSJacob Faibussowitsch PetscCall(MatGolubKahanComputeExplicitOperator(jac->mat[0], jac->B, jac->C, &jac->H, jac->gkbnu)); 11039566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(jac->gkbdelay, &jac->vecz)); 1104e071a0a4SCarola Kruse 1105a51937d4SCarola Kruse ilink = jac->head; 11069566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ilink->ksp, jac->H, jac->H)); 11079566063dSJacob Faibussowitsch if (!jac->suboptionsset) PetscCall(KSPSetFromOptions(ilink->ksp)); 1108e071a0a4SCarola Kruse /* Create gkb_monitor context */ 1109de482cd7SCarola Kruse if (jac->gkbmonitor) { 1110de482cd7SCarola Kruse PetscInt tablevel; 11119566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &jac->gkbviewer)); 11129566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(jac->gkbviewer, PETSCVIEWERASCII)); 11139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetTabLevel((PetscObject)ilink->ksp, &tablevel)); 11149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(jac->gkbviewer, tablevel)); 11159566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->ksp, (PetscObject)ilink->ksp, 1)); 1116de482cd7SCarola Kruse } 11173b224e63SBarry Smith } else { 111868bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 111997bbdb24SBarry Smith i = 0; 11205a9f2f41SSatish Balay ilink = jac->head; 11215a9f2f41SSatish Balay while (ilink) { 11229566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ilink->ksp, jac->mat[i], jac->pmat[i])); 11233b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 11249566063dSJacob Faibussowitsch if (!jac->suboptionsset) PetscCall(KSPSetFromOptions(ilink->ksp)); 112597bbdb24SBarry Smith i++; 11265a9f2f41SSatish Balay ilink = ilink->next; 11270971522cSBarry Smith } 11283b224e63SBarry Smith } 11293b224e63SBarry Smith 11305ddf11f8SNicolas Barnafi /* Set coordinates to the sub PC objects whenever these are set */ 11315ddf11f8SNicolas Barnafi if (jac->coordinates_set) { 11325ddf11f8SNicolas Barnafi PC pc_coords; 11335ddf11f8SNicolas Barnafi if (jac->type == PC_COMPOSITE_SCHUR) { 11345ddf11f8SNicolas Barnafi // Head is first block. 11359566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->head->ksp, &pc_coords)); 11369566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates(pc_coords, jac->head->dim, jac->head->ndofs, jac->head->coords)); 11375ddf11f8SNicolas Barnafi // Second one is Schur block, but its KSP object is in kspschur. 11389566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->kspschur, &pc_coords)); 11399566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates(pc_coords, jac->head->next->dim, jac->head->next->ndofs, jac->head->next->coords)); 11405ddf11f8SNicolas Barnafi } else if (jac->type == PC_COMPOSITE_GKB) { 11419d3446b2SPierre Jolivet PetscCall(PetscInfo(pc, "Warning: Setting coordinates does nothing for the GKB Fieldpslit preconditioner\n")); 11425ddf11f8SNicolas Barnafi } else { 11435ddf11f8SNicolas Barnafi ilink = jac->head; 11445ddf11f8SNicolas Barnafi while (ilink) { 11459566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ilink->ksp, &pc_coords)); 11469566063dSJacob Faibussowitsch PetscCall(PCSetCoordinates(pc_coords, ilink->dim, ilink->ndofs, ilink->coords)); 11475ddf11f8SNicolas Barnafi ilink = ilink->next; 11485ddf11f8SNicolas Barnafi } 11495ddf11f8SNicolas Barnafi } 11505ddf11f8SNicolas Barnafi } 11515ddf11f8SNicolas Barnafi 1152c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 11533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11540971522cSBarry Smith } 11550971522cSBarry Smith 115673716367SStefano Zampini static PetscErrorCode PCSetUpOnBlocks_FieldSplit_Schur(PC pc) 115773716367SStefano Zampini { 115873716367SStefano Zampini PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 115973716367SStefano Zampini PC_FieldSplitLink ilinkA = jac->head; 116073716367SStefano Zampini KSP kspA = ilinkA->ksp, kspUpper = jac->kspupper; 116173716367SStefano Zampini 116273716367SStefano Zampini PetscFunctionBegin; 116373716367SStefano Zampini if (jac->schurfactorization == PC_FIELDSPLIT_SCHUR_FACT_FULL && kspUpper != kspA) { 116473716367SStefano Zampini PetscCall(KSPSetUp(kspUpper)); 116573716367SStefano Zampini PetscCall(KSPSetUpOnBlocks(kspUpper)); 116673716367SStefano Zampini } 116773716367SStefano Zampini PetscCall(KSPSetUp(kspA)); 116873716367SStefano Zampini PetscCall(KSPSetUpOnBlocks(kspA)); 116973716367SStefano Zampini if (jac->schurpre != PC_FIELDSPLIT_SCHUR_PRE_FULL) { 117073716367SStefano Zampini PetscCall(KSPSetUp(jac->kspschur)); 117173716367SStefano Zampini PetscCall(KSPSetUpOnBlocks(jac->kspschur)); 1172ab1f0642SPierre Jolivet } else if (kspUpper == kspA) { 11735becce15SPierre Jolivet Mat A; 1174ab1f0642SPierre Jolivet PetscInt m, M, N; 1175ab1f0642SPierre Jolivet VecType vtype; 1176ab1f0642SPierre Jolivet PetscMemType mtype; 1177ab1f0642SPierre Jolivet PetscScalar *array; 1178ab1f0642SPierre Jolivet 1179ab1f0642SPierre Jolivet PetscCall(MatGetSize(jac->B, &M, &N)); 1180ab1f0642SPierre Jolivet PetscCall(MatGetLocalSize(jac->B, &m, NULL)); 1181ab1f0642SPierre Jolivet PetscCall(MatGetVecType(jac->B, &vtype)); 1182ab1f0642SPierre Jolivet PetscCall(VecGetArrayAndMemType(ilinkA->x, &array, &mtype)); 1183ab1f0642SPierre Jolivet PetscCall(VecRestoreArrayAndMemType(ilinkA->x, &array)); 1184ab1f0642SPierre Jolivet if (PetscMemTypeHost(mtype) || (!PetscDefined(HAVE_CUDA) && !PetscDefined(HAVE_HIP))) PetscCall(PetscMalloc1(m * (N + 1), &array)); 1185ab1f0642SPierre Jolivet #if PetscDefined(HAVE_CUDA) 1186ab1f0642SPierre Jolivet else if (PetscMemTypeCUDA(mtype)) PetscCallCUDA(cudaMalloc((void **)&array, sizeof(PetscScalar) * m * (N + 1))); 1187ab1f0642SPierre Jolivet #endif 1188ab1f0642SPierre Jolivet #if PetscDefined(HAVE_HIP) 1189ab1f0642SPierre Jolivet else if (PetscMemTypeHIP(mtype)) PetscCallHIP(hipMalloc((void **)&array, sizeof(PetscScalar) * m * (N + 1))); 1190ab1f0642SPierre Jolivet #endif 1191ab1f0642SPierre Jolivet PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)jac->schur), vtype, m, PETSC_DECIDE, M, N + 1, -1, array, &A)); // number of columns of the Schur complement plus one 11925becce15SPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)jac->schur, "AinvB", (PetscObject)A)); 11935becce15SPierre Jolivet PetscCall(MatDestroy(&A)); 119473716367SStefano Zampini } 119573716367SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 119673716367SStefano Zampini } 119773716367SStefano Zampini 119873716367SStefano Zampini static PetscErrorCode PCSetUpOnBlocks_FieldSplit(PC pc) 119973716367SStefano Zampini { 120073716367SStefano Zampini PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 120173716367SStefano Zampini PC_FieldSplitLink ilink = jac->head; 120273716367SStefano Zampini 120373716367SStefano Zampini PetscFunctionBegin; 120473716367SStefano Zampini while (ilink) { 120573716367SStefano Zampini PetscCall(KSPSetUp(ilink->ksp)); 120673716367SStefano Zampini PetscCall(KSPSetUpOnBlocks(ilink->ksp)); 120773716367SStefano Zampini ilink = ilink->next; 120873716367SStefano Zampini } 120973716367SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 121073716367SStefano Zampini } 121173716367SStefano Zampini 121273716367SStefano Zampini static PetscErrorCode PCSetUpOnBlocks_FieldSplit_GKB(PC pc) 121373716367SStefano Zampini { 121473716367SStefano Zampini PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 121573716367SStefano Zampini PC_FieldSplitLink ilinkA = jac->head; 121673716367SStefano Zampini KSP ksp = ilinkA->ksp; 121773716367SStefano Zampini 121873716367SStefano Zampini PetscFunctionBegin; 121973716367SStefano Zampini PetscCall(KSPSetUp(ksp)); 122073716367SStefano Zampini PetscCall(KSPSetUpOnBlocks(ksp)); 122173716367SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 122273716367SStefano Zampini } 122373716367SStefano Zampini 1224d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_FieldSplit_Schur(PC pc, Vec x, Vec y) 1225d71ae5a4SJacob Faibussowitsch { 12263b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 12273b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 1228443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 1229d9eadc85SPierre Jolivet Mat AinvB = NULL; 1230ab1f0642SPierre Jolivet PetscInt N, P; 12313b224e63SBarry Smith 12323b224e63SBarry Smith PetscFunctionBegin; 1233c5d2311dSJed Brown switch (jac->schurfactorization) { 1234c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 1235a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 12369566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 12379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 12389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 12399566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 12409566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y)); 12419566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 12429566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 12439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 12459566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1246e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 12479566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y)); 12489566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 1249e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 12509566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 12519566063dSJacob Faibussowitsch PetscCall(VecScale(ilinkD->y, jac->schurscale)); 12529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12539566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 1255c5d2311dSJed Brown break; 1256c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 1257a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 12589566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 12599566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 12609566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 12619566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y)); 12629566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 12639566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 12649566063dSJacob Faibussowitsch PetscCall(MatMult(jac->C, ilinkA->y, ilinkD->x)); 12659566063dSJacob Faibussowitsch PetscCall(VecScale(ilinkD->x, -1.)); 12669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 12679566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12689566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 12699566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1270e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 12719566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y)); 1272e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 12739566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 12749566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 12759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12769566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12779566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 1278c5d2311dSJed Brown break; 1279c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 1280a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 12819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 12829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 12839566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1284e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 12859566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y)); 1286e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 12879566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 12889371c9d4SSatish Balay PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 12899371c9d4SSatish Balay PetscCall(MatMult(jac->B, ilinkD->y, ilinkA->x)); 12909566063dSJacob Faibussowitsch PetscCall(VecScale(ilinkA->x, -1.)); 12919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, ADD_VALUES, SCATTER_FORWARD)); 12929566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12939566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, ADD_VALUES, SCATTER_FORWARD)); 12949566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 12959566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y)); 12969566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 12979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 12989566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 12999566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 13009566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 1301c5d2311dSJed Brown break; 1302c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 1303c238f8cdSStefano Zampini /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1] */ 1304ab1f0642SPierre Jolivet PetscCall(MatGetSize(jac->B, NULL, &P)); 1305ab1f0642SPierre Jolivet N = P; 13069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 13079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 13089566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(KSP_Solve_FS_L, kspLower, ilinkA->x, ilinkA->y, NULL)); 1309d9eadc85SPierre Jolivet if (kspUpper == kspA) { 1310d9eadc85SPierre Jolivet PetscCall(PetscObjectQuery((PetscObject)jac->schur, "AinvB", (PetscObject *)&AinvB)); 1311d9eadc85SPierre Jolivet if (AinvB) { 1312d9eadc85SPierre Jolivet PetscCall(MatGetSize(AinvB, NULL, &N)); 1313ab1f0642SPierre Jolivet if (N > P) { // first time PCApply_FieldSplit_Schur() is called 1314f7cbcdf3SPierre Jolivet PetscMemType mtype; 1315ab1f0642SPierre Jolivet Vec c = NULL; 1316f7cbcdf3SPierre Jolivet PetscScalar *array; 1317ab1f0642SPierre Jolivet PetscInt m, M; 1318d9eadc85SPierre Jolivet 1319ab1f0642SPierre Jolivet PetscCall(MatGetSize(jac->B, &M, NULL)); 1320d9eadc85SPierre Jolivet PetscCall(MatGetLocalSize(jac->B, &m, NULL)); 1321ab1f0642SPierre Jolivet PetscCall(MatDenseGetArrayAndMemType(AinvB, &array, &mtype)); 1322ab1f0642SPierre Jolivet if (PetscMemTypeHost(mtype) || (!PetscDefined(HAVE_CUDA) && !PetscDefined(HAVE_HIP))) PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)jac->schur), 1, m, M, array + m * P, &c)); 1323f7cbcdf3SPierre Jolivet #if PetscDefined(HAVE_CUDA) 1324ab1f0642SPierre Jolivet else if (PetscMemTypeCUDA(mtype)) PetscCall(VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)jac->schur), 1, m, M, array + m * P, &c)); 1325f7cbcdf3SPierre Jolivet #endif 1326f7cbcdf3SPierre Jolivet #if PetscDefined(HAVE_HIP) 1327ab1f0642SPierre Jolivet else if (PetscMemTypeHIP(mtype)) PetscCall(VecCreateMPIHIPWithArray(PetscObjectComm((PetscObject)jac->schur), 1, m, M, array + m * P, &c)); 1328f7cbcdf3SPierre Jolivet #endif 1329ab1f0642SPierre Jolivet PetscCall(MatDenseRestoreArrayAndMemType(AinvB, &array)); 1330ab1f0642SPierre Jolivet PetscCall(VecCopy(ilinkA->x, c)); 1331d9eadc85SPierre Jolivet PetscCall(MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user)); 1332d9eadc85SPierre Jolivet PetscCall(KSPSetOperators(jac->kspschur, jac->schur, jac->schur_user)); 1333f7cbcdf3SPierre Jolivet PetscCall(VecCopy(c, ilinkA->y)); // retrieve the solution as the last column of the composed Mat 1334f7cbcdf3SPierre Jolivet PetscCall(VecDestroy(&c)); 1335d9eadc85SPierre Jolivet } 1336d9eadc85SPierre Jolivet } 1337d9eadc85SPierre Jolivet } 1338ab1f0642SPierre Jolivet if (N == P) PetscCall(KSPSolve(kspLower, ilinkA->x, ilinkA->y)); 13399566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspLower, pc, ilinkA->y)); 13409566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(KSP_Solve_FS_L, kspLower, ilinkA->x, ilinkA->y, NULL)); 13419566063dSJacob Faibussowitsch PetscCall(MatMult(jac->C, ilinkA->y, ilinkD->x)); 13429566063dSJacob Faibussowitsch PetscCall(VecScale(ilinkD->x, -1.0)); 13439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 13449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 13453b224e63SBarry Smith 13469566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1347e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 13489566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->kspschur, ilinkD->x, ilinkD->y)); 1349e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 13509566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 13519566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 13529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 13533b224e63SBarry Smith 1354443836d0SMatthew G Knepley if (kspUpper == kspA) { 1355d9eadc85SPierre Jolivet if (!AinvB) { 13569566063dSJacob Faibussowitsch PetscCall(MatMult(jac->B, ilinkD->y, ilinkA->y)); 13579566063dSJacob Faibussowitsch PetscCall(VecAXPY(ilinkA->x, -1.0, ilinkA->y)); 13589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 13599566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y)); 13609566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 13619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 1362d9eadc85SPierre Jolivet } else PetscCall(MatMultAdd(AinvB, ilinkD->y, ilinkA->y, ilinkA->y)); 1363443836d0SMatthew G Knepley } else { 13649566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 13659566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspA, ilinkA->x, ilinkA->y)); 13669566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 13679566063dSJacob Faibussowitsch PetscCall(MatMult(jac->B, ilinkD->y, ilinkA->x)); 13689566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(KSP_Solve_FS_U, kspUpper, ilinkA->x, ilinkA->z, NULL)); 13699566063dSJacob Faibussowitsch PetscCall(KSPSolve(kspUpper, ilinkA->x, ilinkA->z)); 13709566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(kspUpper, pc, ilinkA->z)); 13719566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(KSP_Solve_FS_U, kspUpper, ilinkA->x, ilinkA->z, NULL)); 13729566063dSJacob Faibussowitsch PetscCall(VecAXPY(ilinkA->y, -1.0, ilinkA->z)); 1373443836d0SMatthew G Knepley } 13749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 13759566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 13769566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 1377c5d2311dSJed Brown } 13783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13793b224e63SBarry Smith } 13803b224e63SBarry Smith 13817b665727SPierre Jolivet static PetscErrorCode PCApplyTranspose_FieldSplit_Schur(PC pc, Vec x, Vec y) 13827b665727SPierre Jolivet { 13837b665727SPierre Jolivet PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 13847b665727SPierre Jolivet PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 13857b665727SPierre Jolivet KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 13867b665727SPierre Jolivet 13877b665727SPierre Jolivet PetscFunctionBegin; 13887b665727SPierre Jolivet switch (jac->schurfactorization) { 13897b665727SPierre Jolivet case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 13907b665727SPierre Jolivet /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 13917b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 13927b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 13937b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 13947b665727SPierre Jolivet PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 13957b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspA, ilinkA->x, ilinkA->y)); 13967b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 13977b665727SPierre Jolivet PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 13987b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 13997b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 14007b665727SPierre Jolivet PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1401e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 14027b665727SPierre Jolivet PetscCall(KSPSolveTranspose(jac->kspschur, ilinkD->x, ilinkD->y)); 1403e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 14047b665727SPierre Jolivet PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 14057b665727SPierre Jolivet PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 14067b665727SPierre Jolivet PetscCall(VecScale(ilinkD->y, jac->schurscale)); 14077b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14087b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14097b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14107b665727SPierre Jolivet break; 14117b665727SPierre Jolivet case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 14127b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 14137b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 14147b665727SPierre Jolivet PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14157b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspA, ilinkA->x, ilinkA->y)); 14167b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 14177b665727SPierre Jolivet PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14187b665727SPierre Jolivet PetscCall(MatMultTranspose(jac->B, ilinkA->y, ilinkD->x)); 14197b665727SPierre Jolivet PetscCall(VecScale(ilinkD->x, -1.)); 14207b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 14217b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14227b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 14237b665727SPierre Jolivet PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1424e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 14257b665727SPierre Jolivet PetscCall(KSPSolveTranspose(jac->kspschur, ilinkD->x, ilinkD->y)); 1426e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 14277b665727SPierre Jolivet PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 14287b665727SPierre Jolivet PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 14297b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14307b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14317b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14327b665727SPierre Jolivet break; 14337b665727SPierre Jolivet case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 14347b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 14357b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 14367b665727SPierre Jolivet PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1437e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 14387b665727SPierre Jolivet PetscCall(KSPSolveTranspose(jac->kspschur, ilinkD->x, ilinkD->y)); 1439e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 14407b665727SPierre Jolivet PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 14417b665727SPierre Jolivet PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 14427b665727SPierre Jolivet PetscCall(MatMultTranspose(jac->C, ilinkD->y, ilinkA->x)); 14437b665727SPierre Jolivet PetscCall(VecScale(ilinkA->x, -1.)); 14447b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, ADD_VALUES, SCATTER_FORWARD)); 14457b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14467b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, ADD_VALUES, SCATTER_FORWARD)); 14477b665727SPierre Jolivet PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14487b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspA, ilinkA->x, ilinkA->y)); 14497b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 14507b665727SPierre Jolivet PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14517b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14527b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14537b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14547b665727SPierre Jolivet break; 14557b665727SPierre Jolivet case PC_FIELDSPLIT_SCHUR_FACT_FULL: 14567b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 14577b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 14587b665727SPierre Jolivet PetscCall(PetscLogEventBegin(KSP_Solve_FS_U, kspUpper, ilinkA->x, ilinkA->y, NULL)); 14597b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspUpper, ilinkA->x, ilinkA->y)); 14607b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspUpper, pc, ilinkA->y)); 14617b665727SPierre Jolivet PetscCall(PetscLogEventEnd(KSP_Solve_FS_U, kspUpper, ilinkA->x, ilinkA->y, NULL)); 14627b665727SPierre Jolivet PetscCall(MatMultTranspose(jac->B, ilinkA->y, ilinkD->x)); 14637b665727SPierre Jolivet PetscCall(VecScale(ilinkD->x, -1.0)); 14647b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 14657b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, ADD_VALUES, SCATTER_FORWARD)); 14667b665727SPierre Jolivet 14677b665727SPierre Jolivet PetscCall(PetscLogEventBegin(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 1468e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, 1)); 14697b665727SPierre Jolivet PetscCall(KSPSolveTranspose(jac->kspschur, ilinkD->x, ilinkD->y)); 1470e0b7e82fSBarry Smith PetscCall(PetscObjectIncrementTabLevel((PetscObject)kspA, (PetscObject)kspA, -1)); 14717b665727SPierre Jolivet PetscCall(KSPCheckSolve(jac->kspschur, pc, ilinkD->y)); 14727b665727SPierre Jolivet PetscCall(PetscLogEventEnd(KSP_Solve_FS_S, jac->kspschur, ilinkD->x, ilinkD->y, NULL)); 14737b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14747b665727SPierre Jolivet 14757b665727SPierre Jolivet if (kspLower == kspA) { 14767b665727SPierre Jolivet PetscCall(MatMultTranspose(jac->C, ilinkD->y, ilinkA->y)); 14777b665727SPierre Jolivet PetscCall(VecAXPY(ilinkA->x, -1.0, ilinkA->y)); 14787b665727SPierre Jolivet PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14797b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspA, ilinkA->x, ilinkA->y)); 14807b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 14817b665727SPierre Jolivet PetscCall(PetscLogEventEnd(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14827b665727SPierre Jolivet } else { 14837b665727SPierre Jolivet PetscCall(PetscLogEventBegin(ilinkA->event, kspA, ilinkA->x, ilinkA->y, NULL)); 14847b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspA, ilinkA->x, ilinkA->y)); 14857b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspA, pc, ilinkA->y)); 14867b665727SPierre Jolivet PetscCall(MatMultTranspose(jac->C, ilinkD->y, ilinkA->x)); 14877b665727SPierre Jolivet PetscCall(PetscLogEventBegin(KSP_Solve_FS_L, kspLower, ilinkA->x, ilinkA->z, NULL)); 14887b665727SPierre Jolivet PetscCall(KSPSolveTranspose(kspLower, ilinkA->x, ilinkA->z)); 14897b665727SPierre Jolivet PetscCall(KSPCheckSolve(kspLower, pc, ilinkA->z)); 14907b665727SPierre Jolivet PetscCall(PetscLogEventEnd(KSP_Solve_FS_L, kspLower, ilinkA->x, ilinkA->z, NULL)); 14917b665727SPierre Jolivet PetscCall(VecAXPY(ilinkA->y, -1.0, ilinkA->z)); 14927b665727SPierre Jolivet } 14937b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14947b665727SPierre Jolivet PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14957b665727SPierre Jolivet PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 14967b665727SPierre Jolivet } 14977b665727SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 14987b665727SPierre Jolivet } 14997b665727SPierre Jolivet 15005becce15SPierre Jolivet #define FieldSplitSplitSolveAdd(ilink, xx, yy) \ 15015becce15SPierre Jolivet ((PetscErrorCode)(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) || \ 15025becce15SPierre Jolivet 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) || \ 15035becce15SPierre Jolivet VecScatterEnd(ilink->sctx, ilink->y, yy, ADD_VALUES, SCATTER_REVERSE))) 15045becce15SPierre Jolivet 1505d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_FieldSplit(PC pc, Vec x, Vec y) 1506d71ae5a4SJacob Faibussowitsch { 15070971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 15085a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 1509939b8a20SBarry Smith PetscInt cnt, bs; 15100971522cSBarry Smith 15110971522cSBarry Smith PetscFunctionBegin; 151279416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 151380670ca5SBarry Smith PetscBool matnest; 151480670ca5SBarry Smith 151580670ca5SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATNEST, &matnest)); 151680670ca5SBarry Smith if (jac->defaultsplit && !matnest) { 15179566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(x, &bs)); 15182472a847SBarry 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); 15199566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(y, &bs)); 15202472a847SBarry 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); 15219566063dSJacob Faibussowitsch PetscCall(VecStrideGatherAll(x, jac->x, INSERT_VALUES)); 15225a9f2f41SSatish Balay while (ilink) { 15239566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15249566063dSJacob Faibussowitsch PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y)); 15259566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y)); 15269566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15275a9f2f41SSatish Balay ilink = ilink->next; 15280971522cSBarry Smith } 15299566063dSJacob Faibussowitsch PetscCall(VecStrideScatterAll(jac->y, y, INSERT_VALUES)); 15301b9fc7fcSBarry Smith } else { 15319566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 15325a9f2f41SSatish Balay while (ilink) { 15339566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAdd(ilink, x, y)); 15345a9f2f41SSatish Balay ilink = ilink->next; 15351b9fc7fcSBarry Smith } 15361b9fc7fcSBarry Smith } 1537e52d2c62SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) { 15389566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 1539e52d2c62SBarry Smith /* solve on first block for first block variables */ 15409566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, INSERT_VALUES, SCATTER_FORWARD)); 15419566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, INSERT_VALUES, SCATTER_FORWARD)); 15429566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15439566063dSJacob Faibussowitsch PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y)); 15449566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y)); 15459566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15469566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 15479566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 1548e52d2c62SBarry Smith 1549e52d2c62SBarry Smith /* compute the residual only onto second block variables using first block variables */ 15509566063dSJacob Faibussowitsch PetscCall(MatMult(jac->Afield[1], ilink->y, ilink->next->x)); 1551e52d2c62SBarry Smith ilink = ilink->next; 15529566063dSJacob Faibussowitsch PetscCall(VecScale(ilink->x, -1.0)); 15539566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD)); 15549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD)); 1555e52d2c62SBarry Smith 1556e52d2c62SBarry Smith /* solve on second block variables */ 15579566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15589566063dSJacob Faibussowitsch PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y)); 15599566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y)); 15609566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15619566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 15629566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 156316913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 156479416396SBarry Smith if (!jac->w1) { 15659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &jac->w1)); 15669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &jac->w2)); 156779416396SBarry Smith } 15689566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 15699566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAdd(ilink, x, y)); 15703e197d65SBarry Smith cnt = 1; 15715a9f2f41SSatish Balay while (ilink->next) { 15725a9f2f41SSatish Balay ilink = ilink->next; 15733e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 15749566063dSJacob Faibussowitsch PetscCall(MatMult(jac->Afield[cnt++], y, ilink->x)); 15759566063dSJacob Faibussowitsch PetscCall(VecScale(ilink->x, -1.0)); 15769566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD)); 15779566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD)); 15789566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15799566063dSJacob Faibussowitsch PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y)); 15809566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y)); 15819566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15829566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 15839566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 15843e197d65SBarry Smith } 158551f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 158611755939SBarry Smith cnt -= 2; 158751f519a2SBarry Smith while (ilink->previous) { 158851f519a2SBarry Smith ilink = ilink->previous; 158911755939SBarry Smith /* compute the residual only over the part of the vector needed */ 15909566063dSJacob Faibussowitsch PetscCall(MatMult(jac->Afield[cnt--], y, ilink->x)); 15919566063dSJacob Faibussowitsch PetscCall(VecScale(ilink->x, -1.0)); 15929566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD)); 15939566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, x, ilink->x, ADD_VALUES, SCATTER_FORWARD)); 15949566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15959566063dSJacob Faibussowitsch PetscCall(KSPSolve(ilink->ksp, ilink->x, ilink->y)); 15969566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y)); 15979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 15989566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 15999566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilink->sctx, ilink->y, y, ADD_VALUES, SCATTER_REVERSE)); 160051f519a2SBarry Smith } 160111755939SBarry Smith } 160263a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unsupported or unknown composition %d", (int)jac->type); 16033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16040971522cSBarry Smith } 16050971522cSBarry Smith 1606d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_FieldSplit_GKB(PC pc, Vec x, Vec y) 1607d71ae5a4SJacob Faibussowitsch { 1608a51937d4SCarola Kruse PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 1609a51937d4SCarola Kruse PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 1610a51937d4SCarola Kruse KSP ksp = ilinkA->ksp; 1611de482cd7SCarola Kruse Vec u, v, Hu, d, work1, work2; 1612e071a0a4SCarola Kruse PetscScalar alpha, z, nrmz2, *vecz; 1613e071a0a4SCarola Kruse PetscReal lowbnd, nu, beta; 1614a51937d4SCarola Kruse PetscInt j, iterGKB; 1615a51937d4SCarola Kruse 1616a51937d4SCarola Kruse PetscFunctionBegin; 16179566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 16189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 16199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, x, ilinkA->x, INSERT_VALUES, SCATTER_FORWARD)); 16209566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, x, ilinkD->x, INSERT_VALUES, SCATTER_FORWARD)); 1621e071a0a4SCarola Kruse 1622e071a0a4SCarola Kruse u = jac->u; 1623e071a0a4SCarola Kruse v = jac->v; 1624e071a0a4SCarola Kruse Hu = jac->Hu; 1625e071a0a4SCarola Kruse d = jac->d; 1626e071a0a4SCarola Kruse work1 = jac->w1; 1627e071a0a4SCarola Kruse work2 = jac->w2; 1628e071a0a4SCarola Kruse vecz = jac->vecz; 1629a51937d4SCarola Kruse 1630a51937d4SCarola Kruse /* Change RHS to comply with matrix regularization H = A + nu*B*B' */ 1631a51937d4SCarola Kruse /* Add q = q + nu*B*b */ 1632a51937d4SCarola Kruse if (jac->gkbnu) { 1633a51937d4SCarola Kruse nu = jac->gkbnu; 16349566063dSJacob Faibussowitsch PetscCall(VecScale(ilinkD->x, jac->gkbnu)); 16359566063dSJacob Faibussowitsch PetscCall(MatMultAdd(jac->B, ilinkD->x, ilinkA->x, ilinkA->x)); /* q = q + nu*B*b */ 1636a51937d4SCarola Kruse } else { 1637a51937d4SCarola Kruse /* Situation when no augmented Lagrangian is used. Then we set inner */ 1638a51937d4SCarola Kruse /* matrix N = I in [Ar13], and thus nu = 1. */ 1639a51937d4SCarola Kruse nu = 1; 1640a51937d4SCarola Kruse } 1641a51937d4SCarola Kruse 1642a51937d4SCarola Kruse /* Transform rhs from [q,tilde{b}] to [0,b] */ 16439566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, ksp, ilinkA->x, ilinkA->y, NULL)); 16449566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, ilinkA->x, ilinkA->y)); 16459566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ksp, pc, ilinkA->y)); 16469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, ksp, ilinkA->x, ilinkA->y, NULL)); 16479566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->B, ilinkA->y, work1)); 16489566063dSJacob Faibussowitsch PetscCall(VecAXPBY(work1, 1.0 / nu, -1.0, ilinkD->x)); /* c = b - B'*x */ 1649a51937d4SCarola Kruse 1650a51937d4SCarola Kruse /* First step of algorithm */ 16519566063dSJacob Faibussowitsch PetscCall(VecNorm(work1, NORM_2, &beta)); /* beta = sqrt(nu*c'*c)*/ 1652e071a0a4SCarola Kruse KSPCheckDot(ksp, beta); 1653addd1e01SJunchao Zhang beta = PetscSqrtReal(nu) * beta; 16549566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v, nu / beta, 0.0, work1)); /* v = nu/beta *c */ 16559566063dSJacob Faibussowitsch PetscCall(MatMult(jac->B, v, work2)); /* u = H^{-1}*B*v */ 16569566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, ksp, work2, u, NULL)); 16579566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, work2, u)); 16589566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ksp, pc, u)); 16599566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, ksp, work2, u, NULL)); 16609566063dSJacob Faibussowitsch PetscCall(MatMult(jac->H, u, Hu)); /* alpha = u'*H*u */ 16619566063dSJacob Faibussowitsch PetscCall(VecDot(Hu, u, &alpha)); 1662e071a0a4SCarola Kruse KSPCheckDot(ksp, alpha); 166308401ef6SPierre Jolivet PetscCheck(PetscRealPart(alpha) > 0.0, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "GKB preconditioner diverged, H is not positive definite"); 1664addd1e01SJunchao Zhang alpha = PetscSqrtReal(PetscAbsScalar(alpha)); 16659566063dSJacob Faibussowitsch PetscCall(VecScale(u, 1.0 / alpha)); 16669566063dSJacob Faibussowitsch PetscCall(VecAXPBY(d, 1.0 / alpha, 0.0, v)); /* v = nu/beta *c */ 1667de482cd7SCarola Kruse 1668a51937d4SCarola Kruse z = beta / alpha; 1669a51937d4SCarola Kruse vecz[1] = z; 1670a51937d4SCarola Kruse 1671de482cd7SCarola Kruse /* Computation of first iterate x(1) and p(1) */ 16729566063dSJacob Faibussowitsch PetscCall(VecAXPY(ilinkA->y, z, u)); 16739566063dSJacob Faibussowitsch PetscCall(VecCopy(d, ilinkD->y)); 16749566063dSJacob Faibussowitsch PetscCall(VecScale(ilinkD->y, -z)); 1675a51937d4SCarola Kruse 16769371c9d4SSatish Balay iterGKB = 1; 16779371c9d4SSatish Balay lowbnd = 2 * jac->gkbtol; 167848a46eb9SPierre Jolivet if (jac->gkbmonitor) PetscCall(PetscViewerASCIIPrintf(jac->gkbviewer, "%3" PetscInt_FMT " GKB Lower bound estimate %14.12e\n", iterGKB, (double)lowbnd)); 1679de482cd7SCarola Kruse 1680a51937d4SCarola Kruse while (iterGKB < jac->gkbmaxit && lowbnd > jac->gkbtol) { 1681a51937d4SCarola Kruse iterGKB += 1; 16829566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->B, u, work1)); /* v <- nu*(B'*u-alpha/nu*v) */ 16839566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v, nu, -alpha, work1)); 16849566063dSJacob Faibussowitsch PetscCall(VecNorm(v, NORM_2, &beta)); /* beta = sqrt(nu)*v'*v */ 1685addd1e01SJunchao Zhang beta = beta / PetscSqrtReal(nu); 16869566063dSJacob Faibussowitsch PetscCall(VecScale(v, 1.0 / beta)); 16879566063dSJacob Faibussowitsch PetscCall(MatMult(jac->B, v, work2)); /* u <- H^{-1}*(B*v-beta*H*u) */ 16889566063dSJacob Faibussowitsch PetscCall(MatMult(jac->H, u, Hu)); 16899566063dSJacob Faibussowitsch PetscCall(VecAXPY(work2, -beta, Hu)); 16909566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilinkA->event, ksp, work2, u, NULL)); 16919566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, work2, u)); 16929566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ksp, pc, u)); 16939566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilinkA->event, ksp, work2, u, NULL)); 16949566063dSJacob Faibussowitsch PetscCall(MatMult(jac->H, u, Hu)); /* alpha = u'*H*u */ 16959566063dSJacob Faibussowitsch PetscCall(VecDot(Hu, u, &alpha)); 1696e071a0a4SCarola Kruse KSPCheckDot(ksp, alpha); 169708401ef6SPierre Jolivet PetscCheck(PetscRealPart(alpha) > 0.0, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "GKB preconditioner diverged, H is not positive definite"); 1698addd1e01SJunchao Zhang alpha = PetscSqrtReal(PetscAbsScalar(alpha)); 16999566063dSJacob Faibussowitsch PetscCall(VecScale(u, 1.0 / alpha)); 1700a51937d4SCarola Kruse 1701e071a0a4SCarola Kruse z = -beta / alpha * z; /* z <- beta/alpha*z */ 1702a51937d4SCarola Kruse vecz[0] = z; 1703a51937d4SCarola Kruse 1704a51937d4SCarola Kruse /* Computation of new iterate x(i+1) and p(i+1) */ 17059566063dSJacob Faibussowitsch PetscCall(VecAXPBY(d, 1.0 / alpha, -beta / alpha, v)); /* d = (v-beta*d)/alpha */ 17069566063dSJacob Faibussowitsch PetscCall(VecAXPY(ilinkA->y, z, u)); /* r = r + z*u */ 17079566063dSJacob Faibussowitsch PetscCall(VecAXPY(ilinkD->y, -z, d)); /* p = p - z*d */ 17089566063dSJacob Faibussowitsch PetscCall(MatMult(jac->H, ilinkA->y, Hu)); /* ||u||_H = u'*H*u */ 17099566063dSJacob Faibussowitsch PetscCall(VecDot(Hu, ilinkA->y, &nrmz2)); 1710a51937d4SCarola Kruse 1711a51937d4SCarola Kruse /* Compute Lower Bound estimate */ 1712a51937d4SCarola Kruse if (iterGKB > jac->gkbdelay) { 1713a51937d4SCarola Kruse lowbnd = 0.0; 1714ad540459SPierre Jolivet for (j = 0; j < jac->gkbdelay; j++) lowbnd += PetscAbsScalar(vecz[j] * vecz[j]); 1715addd1e01SJunchao Zhang lowbnd = PetscSqrtReal(lowbnd / PetscAbsScalar(nrmz2)); 1716a51937d4SCarola Kruse } 1717a51937d4SCarola Kruse 1718ad540459SPierre Jolivet for (j = 0; j < jac->gkbdelay - 1; j++) vecz[jac->gkbdelay - j - 1] = vecz[jac->gkbdelay - j - 2]; 171948a46eb9SPierre Jolivet if (jac->gkbmonitor) PetscCall(PetscViewerASCIIPrintf(jac->gkbviewer, "%3" PetscInt_FMT " GKB Lower bound estimate %14.12e\n", iterGKB, (double)lowbnd)); 1720a51937d4SCarola Kruse } 1721a51937d4SCarola Kruse 17229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 17239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkA->sctx, ilinkA->y, y, INSERT_VALUES, SCATTER_REVERSE)); 17249566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 17259566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ilinkD->sctx, ilinkD->y, y, INSERT_VALUES, SCATTER_REVERSE)); 17263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1727a51937d4SCarola Kruse } 1728a51937d4SCarola Kruse 1729421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink, xx, yy) \ 17303ba16761SJacob Faibussowitsch ((PetscErrorCode)(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) || \ 17319371c9d4SSatish 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) || \ 17323ba16761SJacob Faibussowitsch VecScatterEnd(ilink->sctx, ilink->x, yy, ADD_VALUES, SCATTER_REVERSE))) 1733421e10b8SBarry Smith 1734d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc, Vec x, Vec y) 1735d71ae5a4SJacob Faibussowitsch { 1736421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 1737421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 1738939b8a20SBarry Smith PetscInt bs; 1739421e10b8SBarry Smith 1740421e10b8SBarry Smith PetscFunctionBegin; 1741421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 174280670ca5SBarry Smith PetscBool matnest; 174380670ca5SBarry Smith 174480670ca5SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATNEST, &matnest)); 174580670ca5SBarry Smith if (jac->defaultsplit && !matnest) { 17469566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(x, &bs)); 17472472a847SBarry 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); 17489566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(y, &bs)); 17492472a847SBarry 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); 17509566063dSJacob Faibussowitsch PetscCall(VecStrideGatherAll(x, jac->x, INSERT_VALUES)); 1751421e10b8SBarry Smith while (ilink) { 17529566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 17539566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ilink->ksp, ilink->x, ilink->y)); 17549566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(ilink->ksp, pc, ilink->y)); 17559566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ilink->event, ilink->ksp, ilink->x, ilink->y, NULL)); 1756421e10b8SBarry Smith ilink = ilink->next; 1757421e10b8SBarry Smith } 17589566063dSJacob Faibussowitsch PetscCall(VecStrideScatterAll(jac->y, y, INSERT_VALUES)); 1759421e10b8SBarry Smith } else { 17609566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 1761421e10b8SBarry Smith while (ilink) { 17629566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAddTranspose(ilink, x, y)); 1763421e10b8SBarry Smith ilink = ilink->next; 1764421e10b8SBarry Smith } 1765421e10b8SBarry Smith } 1766421e10b8SBarry Smith } else { 1767421e10b8SBarry Smith if (!jac->w1) { 17689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &jac->w1)); 17699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &jac->w2)); 1770421e10b8SBarry Smith } 17719566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 1772421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 17739566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAddTranspose(ilink, x, y)); 1774421e10b8SBarry Smith while (ilink->next) { 1775421e10b8SBarry Smith ilink = ilink->next; 17769566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, y, jac->w1)); 17779566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->w2, -1.0, jac->w1, x)); 17789566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAddTranspose(ilink, jac->w2, y)); 1779421e10b8SBarry Smith } 1780421e10b8SBarry Smith while (ilink->previous) { 1781421e10b8SBarry Smith ilink = ilink->previous; 17829566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, y, jac->w1)); 17839566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->w2, -1.0, jac->w1, x)); 17849566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAddTranspose(ilink, jac->w2, y)); 1785421e10b8SBarry Smith } 1786421e10b8SBarry Smith } else { 1787421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 1788421e10b8SBarry Smith ilink = ilink->next; 1789421e10b8SBarry Smith } 17909566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAddTranspose(ilink, x, y)); 1791421e10b8SBarry Smith while (ilink->previous) { 1792421e10b8SBarry Smith ilink = ilink->previous; 17939566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, y, jac->w1)); 17949566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->w2, -1.0, jac->w1, x)); 17959566063dSJacob Faibussowitsch PetscCall(FieldSplitSplitSolveAddTranspose(ilink, jac->w2, y)); 1796421e10b8SBarry Smith } 1797421e10b8SBarry Smith } 1798421e10b8SBarry Smith } 17993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1800421e10b8SBarry Smith } 1801421e10b8SBarry Smith 1802d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_FieldSplit(PC pc) 1803d71ae5a4SJacob Faibussowitsch { 18040971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 18055a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head, next; 18060971522cSBarry Smith 18070971522cSBarry Smith PetscFunctionBegin; 18085a9f2f41SSatish Balay while (ilink) { 18099566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ilink->ksp)); 18109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ilink->x)); 18119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ilink->y)); 18129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ilink->z)); 18139566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ilink->sctx)); 18149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilink->is)); 18159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilink->is_col)); 18169566063dSJacob Faibussowitsch PetscCall(PetscFree(ilink->splitname)); 18179566063dSJacob Faibussowitsch PetscCall(PetscFree(ilink->fields)); 18189566063dSJacob Faibussowitsch PetscCall(PetscFree(ilink->fields_col)); 18195a9f2f41SSatish Balay next = ilink->next; 18209566063dSJacob Faibussowitsch PetscCall(PetscFree(ilink)); 18215a9f2f41SSatish Balay ilink = next; 18220971522cSBarry Smith } 1823f5f0d762SBarry Smith jac->head = NULL; 18249566063dSJacob Faibussowitsch PetscCall(PetscFree2(jac->x, jac->y)); 1825574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 18269566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(jac->nsplits, &jac->mat)); 1827574deadeSBarry Smith } else if (jac->mat) { 18280298fd71SBarry Smith jac->mat = NULL; 1829574deadeSBarry Smith } 18309566063dSJacob Faibussowitsch if (jac->pmat) PetscCall(MatDestroyMatrices(jac->nsplits, &jac->pmat)); 18319566063dSJacob Faibussowitsch if (jac->Afield) PetscCall(MatDestroyMatrices(jac->nsplits, &jac->Afield)); 1832f5f0d762SBarry Smith jac->nsplits = 0; 18339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->w1)); 18349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->w2)); 183573716367SStefano Zampini if (jac->schur) PetscCall(PetscObjectCompose((PetscObject)jac->schur, "AinvB", NULL)); 18369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->schur)); 18379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->schurp)); 18389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->schur_user)); 18399566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&jac->kspschur)); 18409566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&jac->kspupper)); 18419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->B)); 18429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->C)); 18439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->H)); 18449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->u)); 18459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->v)); 18469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->Hu)); 18479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->d)); 18489566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->vecz)); 18499566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&jac->gkbviewer)); 18506dbb499eSCian Wilson jac->isrestrict = PETSC_FALSE; 18513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1852574deadeSBarry Smith } 1853574deadeSBarry Smith 1854d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1855d71ae5a4SJacob Faibussowitsch { 1856574deadeSBarry Smith PetscFunctionBegin; 18579566063dSJacob Faibussowitsch PetscCall(PCReset_FieldSplit(pc)); 18589566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 18592e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 18609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetFields_C", NULL)); 18619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetIS_C", NULL)); 18629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetType_C", NULL)); 18639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetBlockSize_C", NULL)); 18642e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitRestrictIS_C", NULL)); 18652e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSchurGetSubKSP_C", NULL)); 18662e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", NULL)); 18672e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBTol_C", NULL)); 18682e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBMaxit_C", NULL)); 18692e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBNu_C", NULL)); 18702e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBDelay_C", NULL)); 18719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurPre_C", NULL)); 18729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSchurPre_C", NULL)); 18739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurFactType_C", NULL)); 18742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurScale_C", NULL)); 18753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18760971522cSBarry Smith } 18770971522cSBarry Smith 1878ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc, PetscOptionItems PetscOptionsObject) 1879d71ae5a4SJacob Faibussowitsch { 18806c924f48SJed Brown PetscInt bs; 18817b752e3dSPatrick Sanan PetscBool flg; 18829dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 18833b224e63SBarry Smith PCCompositeType ctype; 18841b9fc7fcSBarry Smith 18850971522cSBarry Smith PetscFunctionBegin; 1886d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "FieldSplit options"); 18879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_fieldsplit_dm_splits", "Whether to use DMCreateFieldDecomposition() for splits", "PCFieldSplitSetDMSplits", jac->dm_splits, &jac->dm_splits, NULL)); 18889566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_fieldsplit_block_size", "Blocksize that defines number of fields", "PCFieldSplitSetBlockSize", jac->bs, &bs, &flg)); 18891baa6e33SBarry Smith if (flg) PetscCall(PCFieldSplitSetBlockSize(pc, bs)); 18902686e3e9SMatthew G. Knepley jac->diag_use_amat = pc->useAmat; 18919566063dSJacob 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)); 18922686e3e9SMatthew G. Knepley jac->offdiag_use_amat = pc->useAmat; 18939566063dSJacob 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)); 18949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_fieldsplit_detect_saddle_point", "Form 2-way split by detecting zero diagonal entries", "PCFieldSplitSetDetectSaddlePoint", jac->detect, &jac->detect, NULL)); 18959566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetDetectSaddlePoint(pc, jac->detect)); /* Sets split type and Schur PC type */ 18969566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_fieldsplit_type", "Type of composition", "PCFieldSplitSetType", PCCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&ctype, &flg)); 18971baa6e33SBarry Smith if (flg) PetscCall(PCFieldSplitSetType(pc, ctype)); 1898c30613efSMatthew Knepley /* Only setup fields once */ 1899b6555650SPierre Jolivet if (jac->bs > 0 && jac->nsplits == 0) { 190080670ca5SBarry Smith /* only allow user to set fields from command line. 1901d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 19029566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetRuntimeSplits_Private(pc)); 19039566063dSJacob Faibussowitsch if (jac->splitdefined) PetscCall(PetscInfo(pc, "Splits defined using the options database\n")); 1904d32f9abdSBarry Smith } 1905c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 19069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_fieldsplit_schur_factorization_type", PCFieldSplitSchurFactTypes, (PetscEnum *)&jac->schurfactorization, &flg)); 19079566063dSJacob Faibussowitsch if (flg) PetscCall(PetscInfo(pc, "Deprecated use of -pc_fieldsplit_schur_factorization_type\n")); 19089566063dSJacob 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)); 19099566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_fieldsplit_schur_precondition", "How to build preconditioner for Schur complement", "PCFieldSplitSetSchurPre", PCFieldSplitSchurPreTypes, (PetscEnum)jac->schurpre, (PetscEnum *)&jac->schurpre, NULL)); 19109566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-pc_fieldsplit_schur_scale", "Scale Schur complement", "PCFieldSplitSetSchurScale", jac->schurscale, &jac->schurscale, NULL)); 1911a51937d4SCarola Kruse } else if (jac->type == PC_COMPOSITE_GKB) { 1912a077d33dSBarry Smith PetscCall(PetscOptionsReal("-pc_fieldsplit_gkb_tol", "The tolerance for the lower bound stopping criterion", "PCFieldSplitSetGKBTol", jac->gkbtol, &jac->gkbtol, NULL)); 1913a077d33dSBarry Smith PetscCall(PetscOptionsInt("-pc_fieldsplit_gkb_delay", "The delay value for lower bound criterion", "PCFieldSplitSetGKBDelay", jac->gkbdelay, &jac->gkbdelay, NULL)); 1914a077d33dSBarry Smith PetscCall(PetscOptionsBoundedReal("-pc_fieldsplit_gkb_nu", "Parameter in augmented Lagrangian approach", "PCFieldSplitSetGKBNu", jac->gkbnu, &jac->gkbnu, NULL, 0.0)); 1915a077d33dSBarry Smith PetscCall(PetscOptionsInt("-pc_fieldsplit_gkb_maxit", "Maximum allowed number of iterations", "PCFieldSplitSetGKBMaxit", jac->gkbmaxit, &jac->gkbmaxit, NULL)); 19169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_fieldsplit_gkb_monitor", "Prints number of GKB iterations and error", "PCFieldSplitGKB", jac->gkbmonitor, &jac->gkbmonitor, NULL)); 1917c5d2311dSJed Brown } 191834603f55SBarry Smith /* 191934603f55SBarry Smith In the initial call to this routine the sub-solver data structures do not exist so we cannot call KSPSetFromOptions() on them yet. 192034603f55SBarry 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 192134603f55SBarry Smith is called on the outer solver in case changes were made in the options database 192234603f55SBarry Smith 192334603f55SBarry 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() 192434603f55SBarry 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. 192534603f55SBarry Smith Without this extra check test p2p1fetidp_olof_full and others fail with incorrect matrix types. 192634603f55SBarry Smith 192734603f55SBarry Smith There could be a negative side effect of calling the KSPSetFromOptions() below. 192834603f55SBarry Smith 192934603f55SBarry 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 193034603f55SBarry Smith */ 193134603f55SBarry Smith if (jac->issetup) { 193234603f55SBarry Smith PC_FieldSplitLink ilink = jac->head; 193334603f55SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 193434603f55SBarry Smith if (jac->kspupper && jac->kspupper->totalits > 0) PetscCall(KSPSetFromOptions(jac->kspupper)); 193534603f55SBarry Smith if (jac->kspschur && jac->kspschur->totalits > 0) PetscCall(KSPSetFromOptions(jac->kspschur)); 193634603f55SBarry Smith } 193734603f55SBarry Smith while (ilink) { 193834603f55SBarry Smith if (ilink->ksp->totalits > 0) PetscCall(KSPSetFromOptions(ilink->ksp)); 193934603f55SBarry Smith ilink = ilink->next; 194034603f55SBarry Smith } 194134603f55SBarry Smith } 1942d0609cedSBarry Smith PetscOptionsHeadEnd(); 19433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19440971522cSBarry Smith } 19450971522cSBarry Smith 1946d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc, const char splitname[], PetscInt n, const PetscInt *fields, const PetscInt *fields_col) 1947d71ae5a4SJacob Faibussowitsch { 194897bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 19495a9f2f41SSatish Balay PC_FieldSplitLink ilink, next = jac->head; 195069a612a9SBarry Smith char prefix[128]; 19515d4c12cdSJungho Lee PetscInt i; 1952835f2295SStefano Zampini PetscLogEvent nse; 19530971522cSBarry Smith 19540971522cSBarry Smith PetscFunctionBegin; 19556c924f48SJed Brown if (jac->splitdefined) { 19569566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Ignoring new split \"%s\" because the splits have already been defined\n", splitname)); 19573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19586c924f48SJed Brown } 1959ac530a7eSPierre Jolivet for (i = 0; i < n; i++) PetscCheck(fields[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %" PetscInt_FMT " requested", fields[i]); 19609566063dSJacob Faibussowitsch PetscCall(PetscNew(&ilink)); 1961a04f6461SBarry Smith if (splitname) { 19629566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &ilink->splitname)); 1963a04f6461SBarry Smith } else { 19649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3, &ilink->splitname)); 196563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ilink->splitname, 2, "%" PetscInt_FMT, jac->nsplits)); 1966a04f6461SBarry Smith } 1967835f2295SStefano Zampini PetscCall(PetscMPIIntCast(jac->nsplits, &nse)); 1968835f2295SStefano Zampini ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + nse : KSP_Solve_FS_0 + 4; /* Splits greater than 4 logged in 4th split */ 19699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &ilink->fields)); 19709566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ilink->fields, fields, n)); 19719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &ilink->fields_col)); 19729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ilink->fields_col, fields_col, n)); 19732fa5cd67SKarl Rupp 19745a9f2f41SSatish Balay ilink->nfields = n; 19750298fd71SBarry Smith ilink->next = NULL; 19769566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &ilink->ksp)); 19773821be0aSBarry Smith PetscCall(KSPSetNestLevel(ilink->ksp, pc->kspnestlevel)); 19789566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ilink->ksp, pc->erroriffailure)); 19799566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->ksp, (PetscObject)pc, 1)); 19809566063dSJacob Faibussowitsch PetscCall(KSPSetType(ilink->ksp, KSPPREONLY)); 198169a612a9SBarry Smith 19829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname)); 19839566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ilink->ksp, prefix)); 19840971522cSBarry Smith 19850971522cSBarry Smith if (!next) { 19865a9f2f41SSatish Balay jac->head = ilink; 19870298fd71SBarry Smith ilink->previous = NULL; 19880971522cSBarry Smith } else { 1989ad540459SPierre Jolivet while (next->next) next = next->next; 19905a9f2f41SSatish Balay next->next = ilink; 199151f519a2SBarry Smith ilink->previous = next; 19920971522cSBarry Smith } 19930971522cSBarry Smith jac->nsplits++; 19943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19950971522cSBarry Smith } 19960971522cSBarry Smith 1997d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSchurGetSubKSP_FieldSplit(PC pc, PetscInt *n, KSP **subksp) 1998d71ae5a4SJacob Faibussowitsch { 1999285fb4e2SStefano Zampini PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2000285fb4e2SStefano Zampini 2001285fb4e2SStefano Zampini PetscFunctionBegin; 2002285fb4e2SStefano Zampini *subksp = NULL; 2003285fb4e2SStefano Zampini if (n) *n = 0; 2004285fb4e2SStefano Zampini if (jac->type == PC_COMPOSITE_SCHUR) { 2005285fb4e2SStefano Zampini PetscInt nn; 2006285fb4e2SStefano Zampini 200728b400f6SJacob Faibussowitsch PetscCheck(jac->schur, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitSchurGetSubKSP()"); 200863a3b9bcSJacob Faibussowitsch PetscCheck(jac->nsplits == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Unexpected number of splits %" PetscInt_FMT " != 2", jac->nsplits); 2009285fb4e2SStefano Zampini nn = jac->nsplits + (jac->kspupper != jac->head->ksp ? 1 : 0); 20109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nn, subksp)); 2011285fb4e2SStefano Zampini (*subksp)[0] = jac->head->ksp; 2012285fb4e2SStefano Zampini (*subksp)[1] = jac->kspschur; 2013285fb4e2SStefano Zampini if (jac->kspupper != jac->head->ksp) (*subksp)[2] = jac->kspupper; 2014285fb4e2SStefano Zampini if (n) *n = nn; 2015285fb4e2SStefano Zampini } 20163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2017285fb4e2SStefano Zampini } 2018285fb4e2SStefano Zampini 2019d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc, PetscInt *n, KSP **subksp) 2020d71ae5a4SJacob Faibussowitsch { 2021e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2022e69d4d44SBarry Smith 2023e69d4d44SBarry Smith PetscFunctionBegin; 202428b400f6SJacob Faibussowitsch PetscCheck(jac->schur, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()"); 20259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->nsplits, subksp)); 20269566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(jac->schur, *subksp)); 20272fa5cd67SKarl Rupp 2028e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 202913e0d083SBarry Smith if (n) *n = jac->nsplits; 20303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2031e69d4d44SBarry Smith } 20320971522cSBarry Smith 2033d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc, PetscInt *n, KSP **subksp) 2034d71ae5a4SJacob Faibussowitsch { 20350971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 20360971522cSBarry Smith PetscInt cnt = 0; 20375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 20380971522cSBarry Smith 20390971522cSBarry Smith PetscFunctionBegin; 20409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->nsplits, subksp)); 20415a9f2f41SSatish Balay while (ilink) { 20425a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 20435a9f2f41SSatish Balay ilink = ilink->next; 20440971522cSBarry Smith } 204563a3b9bcSJacob 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); 204613e0d083SBarry Smith if (n) *n = jac->nsplits; 20473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20480971522cSBarry Smith } 20490971522cSBarry Smith 2050cc4c1da9SBarry Smith /*@ 2051f1580f4eSBarry Smith PCFieldSplitRestrictIS - Restricts the fieldsplit `IS`s to be within a given `IS`. 20526dbb499eSCian Wilson 20536dbb499eSCian Wilson Input Parameters: 20546dbb499eSCian Wilson + pc - the preconditioner context 2055feefa0e1SJacob Faibussowitsch - isy - the index set that defines the indices to which the fieldsplit is to be restricted 20566dbb499eSCian Wilson 20576dbb499eSCian Wilson Level: advanced 20586dbb499eSCian Wilson 2059feefa0e1SJacob Faibussowitsch Developer Notes: 2060f1580f4eSBarry Smith It seems the resulting `IS`s will not cover the entire space, so 2061f1580f4eSBarry Smith how can they define a convergent preconditioner? Needs explaining. 2062f1580f4eSBarry Smith 206360f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()` 20646dbb499eSCian Wilson @*/ 2065d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitRestrictIS(PC pc, IS isy) 2066d71ae5a4SJacob Faibussowitsch { 20676dbb499eSCian Wilson PetscFunctionBegin; 20686dbb499eSCian Wilson PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20696dbb499eSCian Wilson PetscValidHeaderSpecific(isy, IS_CLASSID, 2); 2070cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitRestrictIS_C", (PC, IS), (pc, isy)); 20713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20726dbb499eSCian Wilson } 20736dbb499eSCian Wilson 2074d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy) 2075d71ae5a4SJacob Faibussowitsch { 20766dbb499eSCian Wilson PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 20776dbb499eSCian Wilson PC_FieldSplitLink ilink = jac->head, next; 20786dbb499eSCian Wilson PetscInt localsize, size, sizez, i; 20796dbb499eSCian Wilson const PetscInt *ind, *indz; 20806dbb499eSCian Wilson PetscInt *indc, *indcz; 20816dbb499eSCian Wilson PetscBool flg; 20826dbb499eSCian Wilson 20836dbb499eSCian Wilson PetscFunctionBegin; 20849566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isy, &localsize)); 20859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&localsize, &size, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)isy))); 20866dbb499eSCian Wilson size -= localsize; 20876dbb499eSCian Wilson while (ilink) { 20886dbb499eSCian Wilson IS isrl, isr; 20891c7cfba8SBarry Smith PC subpc; 20909566063dSJacob Faibussowitsch PetscCall(ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl)); 20919566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrl, &localsize)); 20929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(localsize, &indc)); 20939566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrl, &ind)); 20949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indc, ind, localsize)); 20959566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrl, &ind)); 20969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrl)); 20976dbb499eSCian Wilson for (i = 0; i < localsize; i++) *(indc + i) += size; 20989566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)isy), localsize, indc, PETSC_OWN_POINTER, &isr)); 20999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isr)); 21009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilink->is)); 21016dbb499eSCian Wilson ilink->is = isr; 21029566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isr)); 21039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilink->is_col)); 21046dbb499eSCian Wilson ilink->is_col = isr; 21059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isr)); 21069566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ilink->ksp, &subpc)); 21079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)subpc, PCFIELDSPLIT, &flg)); 21086dbb499eSCian Wilson if (flg) { 21096dbb499eSCian Wilson IS iszl, isz; 21106dbb499eSCian Wilson MPI_Comm comm; 21119566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ilink->is, &localsize)); 21126dbb499eSCian Wilson comm = PetscObjectComm((PetscObject)ilink->is); 21139566063dSJacob Faibussowitsch PetscCall(ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl)); 21149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&localsize, &sizez, 1, MPIU_INT, MPI_SUM, comm)); 21156dbb499eSCian Wilson sizez -= localsize; 21169566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iszl, &localsize)); 21179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(localsize, &indcz)); 21189566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iszl, &indz)); 21199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indcz, indz, localsize)); 21209566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iszl, &indz)); 21219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iszl)); 21226dbb499eSCian Wilson for (i = 0; i < localsize; i++) *(indcz + i) += sizez; 21239566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, localsize, indcz, PETSC_OWN_POINTER, &isz)); 21249566063dSJacob Faibussowitsch PetscCall(PCFieldSplitRestrictIS(subpc, isz)); 21259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isz)); 21266dbb499eSCian Wilson } 21276dbb499eSCian Wilson next = ilink->next; 21286dbb499eSCian Wilson ilink = next; 21296dbb499eSCian Wilson } 21306dbb499eSCian Wilson jac->isrestrict = PETSC_TRUE; 21313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21326dbb499eSCian Wilson } 21336dbb499eSCian Wilson 2134d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc, const char splitname[], IS is) 2135d71ae5a4SJacob Faibussowitsch { 2136704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2137704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 2138704ba839SBarry Smith char prefix[128]; 2139835f2295SStefano Zampini PetscLogEvent nse; 2140704ba839SBarry Smith 2141704ba839SBarry Smith PetscFunctionBegin; 21426c924f48SJed Brown if (jac->splitdefined) { 21439566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Ignoring new split \"%s\" because the splits have already been defined\n", splitname)); 21443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21456c924f48SJed Brown } 21469566063dSJacob Faibussowitsch PetscCall(PetscNew(&ilink)); 2147a04f6461SBarry Smith if (splitname) { 21489566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &ilink->splitname)); 2149a04f6461SBarry Smith } else { 21509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(8, &ilink->splitname)); 215163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(ilink->splitname, 7, "%" PetscInt_FMT, jac->nsplits)); 2152a04f6461SBarry Smith } 2153835f2295SStefano Zampini PetscCall(PetscMPIIntCast(jac->nsplits, &nse)); 2154835f2295SStefano Zampini ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + nse : KSP_Solve_FS_0 + 4; /* Splits greater than 4 logged in 4th split */ 21559566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is)); 21569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilink->is)); 2157b5787286SJed Brown ilink->is = is; 21589566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is)); 21599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilink->is_col)); 2160b5787286SJed Brown ilink->is_col = is; 21610298fd71SBarry Smith ilink->next = NULL; 21629566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &ilink->ksp)); 21633821be0aSBarry Smith PetscCall(KSPSetNestLevel(ilink->ksp, pc->kspnestlevel)); 21649566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ilink->ksp, pc->erroriffailure)); 21659566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->ksp, (PetscObject)pc, 1)); 21669566063dSJacob Faibussowitsch PetscCall(KSPSetType(ilink->ksp, KSPPREONLY)); 2167704ba839SBarry Smith 21689566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname)); 21699566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ilink->ksp, prefix)); 2170704ba839SBarry Smith 2171704ba839SBarry Smith if (!next) { 2172704ba839SBarry Smith jac->head = ilink; 21730298fd71SBarry Smith ilink->previous = NULL; 2174704ba839SBarry Smith } else { 2175ad540459SPierre Jolivet while (next->next) next = next->next; 2176704ba839SBarry Smith next->next = ilink; 2177704ba839SBarry Smith ilink->previous = next; 2178704ba839SBarry Smith } 2179704ba839SBarry Smith jac->nsplits++; 21803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2181704ba839SBarry Smith } 2182704ba839SBarry Smith 21835d83a8b1SBarry Smith /*@ 218460f59c3bSBarry Smith PCFieldSplitSetFields - Sets the fields that define one particular split in `PCFIELDSPLIT` 21850971522cSBarry Smith 2186c3339decSBarry Smith Logically Collective 21870971522cSBarry Smith 21880971522cSBarry Smith Input Parameters: 21890971522cSBarry Smith + pc - the preconditioner context 219060f59c3bSBarry Smith . splitname - name of this split, if `NULL` the number of the split is used 21910971522cSBarry Smith . n - the number of fields in this split 21923b374dbdSBarry Smith . fields - the fields in this split 219380670ca5SBarry Smith - fields_col - generally the same as `fields`, if it does not match `fields` then the submatrix that is solved for this set of fields comes from an off-diagonal block 219480670ca5SBarry Smith of the matrix and `fields_col` provides the column indices for that block 219580670ca5SBarry Smith 219680670ca5SBarry Smith Options Database Key: 219780670ca5SBarry Smith . -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the `%d`'th split 21980971522cSBarry Smith 21990971522cSBarry Smith Level: intermediate 22000971522cSBarry Smith 220195452b02SPatrick Sanan Notes: 2202f1580f4eSBarry Smith Use `PCFieldSplitSetIS()` to set a general set of indices as a split. 2203d32f9abdSBarry Smith 220480670ca5SBarry Smith If the matrix used to construct the preconditioner is `MATNEST` then field i refers to the `is_row[i]` `IS` passed to `MatCreateNest()`. 220580670ca5SBarry Smith 220680670ca5SBarry Smith If the matrix used to construct the preconditioner is not `MATNEST` then 220754c05997SPierre Jolivet `PCFieldSplitSetFields()` is for defining fields as strided blocks (based on the block size provided to the matrix with `MatSetBlockSize()` or 220880670ca5SBarry Smith to the `PC` with `PCFieldSplitSetBlockSize()`). For example, if the block 2209f1580f4eSBarry 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 221014c74629SNuno Nobre 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x23x56x8.. x12x45x78x.... 2211f1580f4eSBarry Smith where the numbered entries indicate what is in the split. 2212d32f9abdSBarry Smith 2213db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 221460f59c3bSBarry Smith for this split will be available under the prefix `-fieldsplit_SPLITNAME_`. 2215db4c96c1SJed Brown 221680670ca5SBarry Smith `PCFieldSplitSetIS()` does not support having a `fields_col` different from `fields` 22173b374dbdSBarry Smith 2218feefa0e1SJacob Faibussowitsch Developer Notes: 2219f1580f4eSBarry Smith This routine does not actually create the `IS` representing the split, that is delayed 2220f1580f4eSBarry Smith until `PCSetUp_FieldSplit()`, because information about the vector/matrix layouts may not be 22215d4c12cdSJungho Lee available when this routine is called. 22225d4c12cdSJungho Lee 222380670ca5SBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetBlockSize()`, `PCFieldSplitSetIS()`, `PCFieldSplitRestrictIS()`, 222454c05997SPierre Jolivet `MatSetBlockSize()`, `MatCreateNest()` 22250971522cSBarry Smith @*/ 2226cc4c1da9SBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc, const char splitname[], PetscInt n, const PetscInt fields[], const PetscInt fields_col[]) 2227d71ae5a4SJacob Faibussowitsch { 22280971522cSBarry Smith PetscFunctionBegin; 22290700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 22304f572ea9SToby Isaac PetscAssertPointer(splitname, 2); 223163a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Provided number of fields %" PetscInt_FMT " in split \"%s\" not positive", n, splitname); 22324f572ea9SToby Isaac PetscAssertPointer(fields, 4); 2233cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetFields_C", (PC, const char[], PetscInt, const PetscInt *, const PetscInt *), (pc, splitname, n, fields, fields_col)); 22343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22350971522cSBarry Smith } 22360971522cSBarry Smith 2237c84da90fSDmitry Karpeev /*@ 2238f1580f4eSBarry Smith PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat) to build 2239c14e9f00SDavid Andrs the sub-matrices associated with each split. Where `KSPSetOperators`(ksp,Amat,Pmat) was used to supply the operators. 2240c84da90fSDmitry Karpeev 2241c3339decSBarry Smith Logically Collective 2242c84da90fSDmitry Karpeev 2243c84da90fSDmitry Karpeev Input Parameters: 2244c84da90fSDmitry Karpeev + pc - the preconditioner object 2245c84da90fSDmitry Karpeev - flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from 2246c84da90fSDmitry Karpeev 224720f4b53cSBarry Smith Options Database Key: 2248147403d9SBarry Smith . -pc_fieldsplit_diag_use_amat - use the Amat to provide the diagonal blocks 2249c84da90fSDmitry Karpeev 2250c84da90fSDmitry Karpeev Level: intermediate 2251c84da90fSDmitry Karpeev 225260f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitGetDiagUseAmat()`, `PCFieldSplitSetOffDiagUseAmat()`, `PCFIELDSPLIT` 2253c84da90fSDmitry Karpeev @*/ 2254d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetDiagUseAmat(PC pc, PetscBool flg) 2255d71ae5a4SJacob Faibussowitsch { 2256c84da90fSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2257c84da90fSDmitry Karpeev PetscBool isfs; 2258c84da90fSDmitry Karpeev 2259c84da90fSDmitry Karpeev PetscFunctionBegin; 2260c84da90fSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 22619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs)); 226228b400f6SJacob Faibussowitsch PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT); 2263c84da90fSDmitry Karpeev jac->diag_use_amat = flg; 22643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2265c84da90fSDmitry Karpeev } 2266c84da90fSDmitry Karpeev 2267c84da90fSDmitry Karpeev /*@ 2268f1580f4eSBarry Smith PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat) to build 2269c14e9f00SDavid Andrs the sub-matrices associated with each split. Where `KSPSetOperators`(ksp,Amat,Pmat) was used to supply the operators. 2270c84da90fSDmitry Karpeev 2271c3339decSBarry Smith Logically Collective 2272c84da90fSDmitry Karpeev 227320f4b53cSBarry Smith Input Parameter: 2274c84da90fSDmitry Karpeev . pc - the preconditioner object 2275c84da90fSDmitry Karpeev 227620f4b53cSBarry Smith Output Parameter: 2277c84da90fSDmitry Karpeev . flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from 2278c84da90fSDmitry Karpeev 2279c84da90fSDmitry Karpeev Level: intermediate 2280c84da90fSDmitry Karpeev 228160f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitSetDiagUseAmat()`, `PCFieldSplitGetOffDiagUseAmat()`, `PCFIELDSPLIT` 2282c84da90fSDmitry Karpeev @*/ 2283d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetDiagUseAmat(PC pc, PetscBool *flg) 2284d71ae5a4SJacob Faibussowitsch { 2285c84da90fSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2286c84da90fSDmitry Karpeev PetscBool isfs; 2287c84da90fSDmitry Karpeev 2288c84da90fSDmitry Karpeev PetscFunctionBegin; 2289c84da90fSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 22904f572ea9SToby Isaac PetscAssertPointer(flg, 2); 22919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs)); 229228b400f6SJacob Faibussowitsch PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT); 2293c84da90fSDmitry Karpeev *flg = jac->diag_use_amat; 22943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2295c84da90fSDmitry Karpeev } 2296c84da90fSDmitry Karpeev 2297c84da90fSDmitry Karpeev /*@ 2298f1580f4eSBarry Smith PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat) to build 2299c14e9f00SDavid Andrs the sub-matrices associated with each split. Where `KSPSetOperators`(ksp,Amat,Pmat) was used to supply the operators. 2300c84da90fSDmitry Karpeev 2301c3339decSBarry Smith Logically Collective 2302c84da90fSDmitry Karpeev 2303c84da90fSDmitry Karpeev Input Parameters: 2304c84da90fSDmitry Karpeev + pc - the preconditioner object 2305c84da90fSDmitry Karpeev - flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from 2306c84da90fSDmitry Karpeev 230720f4b53cSBarry Smith Options Database Key: 2308147403d9SBarry Smith . -pc_fieldsplit_off_diag_use_amat <bool> - use the Amat to extract the off-diagonal blocks 2309c84da90fSDmitry Karpeev 2310c84da90fSDmitry Karpeev Level: intermediate 2311c84da90fSDmitry Karpeev 231260f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitGetOffDiagUseAmat()`, `PCFieldSplitSetDiagUseAmat()`, `PCFIELDSPLIT` 2313c84da90fSDmitry Karpeev @*/ 2314d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC pc, PetscBool flg) 2315d71ae5a4SJacob Faibussowitsch { 2316c84da90fSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2317c84da90fSDmitry Karpeev PetscBool isfs; 2318c84da90fSDmitry Karpeev 2319c84da90fSDmitry Karpeev PetscFunctionBegin; 2320c84da90fSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 23219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs)); 232228b400f6SJacob Faibussowitsch PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT); 2323c84da90fSDmitry Karpeev jac->offdiag_use_amat = flg; 23243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2325c84da90fSDmitry Karpeev } 2326c84da90fSDmitry Karpeev 2327c84da90fSDmitry Karpeev /*@ 2328f1580f4eSBarry Smith PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat) to build 2329c14e9f00SDavid Andrs the sub-matrices associated with each split. Where `KSPSetOperators`(ksp,Amat,Pmat) was used to supply the operators. 2330c84da90fSDmitry Karpeev 2331c3339decSBarry Smith Logically Collective 2332c84da90fSDmitry Karpeev 233320f4b53cSBarry Smith Input Parameter: 2334c84da90fSDmitry Karpeev . pc - the preconditioner object 2335c84da90fSDmitry Karpeev 233620f4b53cSBarry Smith Output Parameter: 2337c84da90fSDmitry Karpeev . flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from 2338c84da90fSDmitry Karpeev 2339c84da90fSDmitry Karpeev Level: intermediate 2340c84da90fSDmitry Karpeev 234160f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCSetOperators()`, `KSPSetOperators()`, `PCFieldSplitSetOffDiagUseAmat()`, `PCFieldSplitGetDiagUseAmat()`, `PCFIELDSPLIT` 2342c84da90fSDmitry Karpeev @*/ 2343d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC pc, PetscBool *flg) 2344d71ae5a4SJacob Faibussowitsch { 2345c84da90fSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2346c84da90fSDmitry Karpeev PetscBool isfs; 2347c84da90fSDmitry Karpeev 2348c84da90fSDmitry Karpeev PetscFunctionBegin; 2349c84da90fSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 23504f572ea9SToby Isaac PetscAssertPointer(flg, 2); 23519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs)); 235228b400f6SJacob Faibussowitsch PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PC not of type %s", PCFIELDSPLIT); 2353c84da90fSDmitry Karpeev *flg = jac->offdiag_use_amat; 23543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2355c84da90fSDmitry Karpeev } 2356c84da90fSDmitry Karpeev 2357cc4c1da9SBarry Smith /*@ 2358f1580f4eSBarry Smith PCFieldSplitSetIS - Sets the exact elements for a split in a `PCFIELDSPLIT` 2359704ba839SBarry Smith 2360c3339decSBarry Smith Logically Collective 2361704ba839SBarry Smith 2362704ba839SBarry Smith Input Parameters: 2363704ba839SBarry Smith + pc - the preconditioner context 236460f59c3bSBarry Smith . splitname - name of this split, if `NULL` the number of the split is used 2365f1580f4eSBarry Smith - is - the index set that defines the elements in this split 2366704ba839SBarry Smith 236760f59c3bSBarry Smith Level: intermediate 236860f59c3bSBarry Smith 2369a6ffb8dbSJed Brown Notes: 237080670ca5SBarry Smith Use `PCFieldSplitSetFields()`, for splits defined by strided `IS` based on the matrix block size or the `is_rows[]` passed into `MATNEST` 2371a6ffb8dbSJed Brown 2372db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 2373db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 2374d32f9abdSBarry Smith 237580670ca5SBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetBlockSize()`, `PCFieldSplitSetFields()` 2376704ba839SBarry Smith @*/ 2377d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetIS(PC pc, const char splitname[], IS is) 2378d71ae5a4SJacob Faibussowitsch { 2379704ba839SBarry Smith PetscFunctionBegin; 23800700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 23814f572ea9SToby Isaac if (splitname) PetscAssertPointer(splitname, 2); 2382db4c96c1SJed Brown PetscValidHeaderSpecific(is, IS_CLASSID, 3); 2383cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetIS_C", (PC, const char[], IS), (pc, splitname, is)); 23843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2385704ba839SBarry Smith } 2386704ba839SBarry Smith 2387cc4c1da9SBarry Smith /*@ 2388f1580f4eSBarry Smith PCFieldSplitGetIS - Retrieves the elements for a split as an `IS` 238957a9adfeSMatthew G Knepley 2390c3339decSBarry Smith Logically Collective 239157a9adfeSMatthew G Knepley 239257a9adfeSMatthew G Knepley Input Parameters: 239357a9adfeSMatthew G Knepley + pc - the preconditioner context 239457a9adfeSMatthew G Knepley - splitname - name of this split 239557a9adfeSMatthew G Knepley 239657a9adfeSMatthew G Knepley Output Parameter: 2397feefa0e1SJacob Faibussowitsch . is - the index set that defines the elements in this split, or `NULL` if the split is not found 239857a9adfeSMatthew G Knepley 239957a9adfeSMatthew G Knepley Level: intermediate 240057a9adfeSMatthew G Knepley 240180670ca5SBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetIS()`, `PCFieldSplitGetISByIndex()` 240257a9adfeSMatthew G Knepley @*/ 2403d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetIS(PC pc, const char splitname[], IS *is) 2404d71ae5a4SJacob Faibussowitsch { 240557a9adfeSMatthew G Knepley PetscFunctionBegin; 240657a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 24074f572ea9SToby Isaac PetscAssertPointer(splitname, 2); 24084f572ea9SToby Isaac PetscAssertPointer(is, 3); 240957a9adfeSMatthew G Knepley { 241057a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 241157a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 241257a9adfeSMatthew G Knepley PetscBool found; 241357a9adfeSMatthew G Knepley 24140298fd71SBarry Smith *is = NULL; 241557a9adfeSMatthew G Knepley while (ilink) { 24169566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ilink->splitname, splitname, &found)); 241757a9adfeSMatthew G Knepley if (found) { 241857a9adfeSMatthew G Knepley *is = ilink->is; 241957a9adfeSMatthew G Knepley break; 242057a9adfeSMatthew G Knepley } 242157a9adfeSMatthew G Knepley ilink = ilink->next; 242257a9adfeSMatthew G Knepley } 242357a9adfeSMatthew G Knepley } 24243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 242557a9adfeSMatthew G Knepley } 242657a9adfeSMatthew G Knepley 2427cc4c1da9SBarry Smith /*@ 2428f1580f4eSBarry Smith PCFieldSplitGetISByIndex - Retrieves the elements for a given split as an `IS` 2429998f007dSPierre Jolivet 2430c3339decSBarry Smith Logically Collective 2431998f007dSPierre Jolivet 2432998f007dSPierre Jolivet Input Parameters: 2433998f007dSPierre Jolivet + pc - the preconditioner context 2434998f007dSPierre Jolivet - index - index of this split 2435998f007dSPierre Jolivet 2436998f007dSPierre Jolivet Output Parameter: 2437feefa0e1SJacob Faibussowitsch . is - the index set that defines the elements in this split 2438998f007dSPierre Jolivet 2439998f007dSPierre Jolivet Level: intermediate 2440998f007dSPierre Jolivet 244180670ca5SBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitGetIS()`, `PCFieldSplitSetIS()`, 244280670ca5SBarry Smith 2443998f007dSPierre Jolivet @*/ 2444d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetISByIndex(PC pc, PetscInt index, IS *is) 2445d71ae5a4SJacob Faibussowitsch { 2446998f007dSPierre Jolivet PetscFunctionBegin; 244763a3b9bcSJacob Faibussowitsch PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %" PetscInt_FMT " requested", index); 2448998f007dSPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 24494f572ea9SToby Isaac PetscAssertPointer(is, 3); 2450998f007dSPierre Jolivet { 2451998f007dSPierre Jolivet PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2452998f007dSPierre Jolivet PC_FieldSplitLink ilink = jac->head; 2453998f007dSPierre Jolivet PetscInt i = 0; 245463a3b9bcSJacob Faibussowitsch PetscCheck(index < jac->nsplits, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " requested but only %" PetscInt_FMT " exist", index, jac->nsplits); 2455998f007dSPierre Jolivet 2456998f007dSPierre Jolivet while (i < index) { 2457998f007dSPierre Jolivet ilink = ilink->next; 2458998f007dSPierre Jolivet ++i; 2459998f007dSPierre Jolivet } 24609566063dSJacob Faibussowitsch PetscCall(PCFieldSplitGetIS(pc, ilink->splitname, is)); 2461998f007dSPierre Jolivet } 24623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2463998f007dSPierre Jolivet } 2464998f007dSPierre Jolivet 246551f519a2SBarry Smith /*@ 246651f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 246780670ca5SBarry Smith fieldsplit preconditioner when calling `PCFieldSplitSetFields()`. If not set the matrix block size is used. 246851f519a2SBarry Smith 2469c3339decSBarry Smith Logically Collective 247051f519a2SBarry Smith 247151f519a2SBarry Smith Input Parameters: 247251f519a2SBarry Smith + pc - the preconditioner context 247351f519a2SBarry Smith - bs - the block size 247451f519a2SBarry Smith 247551f519a2SBarry Smith Level: intermediate 247651f519a2SBarry Smith 247780670ca5SBarry Smith Note: 247880670ca5SBarry Smith If the matrix is a `MATNEST` then the `is_rows[]` passed to `MatCreateNest()` determines the fields. 247980670ca5SBarry Smith 248060f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()` 248151f519a2SBarry Smith @*/ 2482d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetBlockSize(PC pc, PetscInt bs) 2483d71ae5a4SJacob Faibussowitsch { 248451f519a2SBarry Smith PetscFunctionBegin; 24850700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2486c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, bs, 2); 2487cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetBlockSize_C", (PC, PetscInt), (pc, bs)); 24883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 248951f519a2SBarry Smith } 249051f519a2SBarry Smith 24910971522cSBarry Smith /*@C 2492f1580f4eSBarry Smith PCFieldSplitGetSubKSP - Gets the `KSP` contexts for all splits 24930971522cSBarry Smith 2494c3339decSBarry Smith Collective 24950971522cSBarry Smith 24960971522cSBarry Smith Input Parameter: 24970971522cSBarry Smith . pc - the preconditioner context 24980971522cSBarry Smith 24990971522cSBarry Smith Output Parameters: 250013e0d083SBarry Smith + n - the number of splits 2501f1580f4eSBarry Smith - subksp - the array of `KSP` contexts 2502196cc216SBarry Smith 25030971522cSBarry Smith Level: advanced 25040971522cSBarry Smith 2505f1580f4eSBarry Smith Notes: 2506f1580f4eSBarry Smith After `PCFieldSplitGetSubKSP()` the array of `KSP`s is to be freed by the user with `PetscFree()` 2507f1580f4eSBarry Smith (not the `KSP`, just the array that contains them). 2508f1580f4eSBarry Smith 2509f1580f4eSBarry Smith You must call `PCSetUp()` before calling `PCFieldSplitGetSubKSP()`. 2510f1580f4eSBarry Smith 2511f1580f4eSBarry Smith If the fieldsplit is of type `PC_COMPOSITE_SCHUR`, it returns the `KSP` object used inside the 2512f1580f4eSBarry Smith Schur complement and the `KSP` object used to iterate over the Schur complement. 2513f1580f4eSBarry Smith To access all the `KSP` objects used in `PC_COMPOSITE_SCHUR`, use `PCFieldSplitSchurGetSubKSP()`. 2514f1580f4eSBarry Smith 2515f1580f4eSBarry Smith If the fieldsplit is of type `PC_COMPOSITE_GKB`, it returns the `KSP` object used to solve the 2516f1580f4eSBarry Smith inner linear system defined by the matrix H in each loop. 2517f1580f4eSBarry Smith 2518feaf08eaSBarry Smith Fortran Note: 2519e41f517fSBarry Smith Call `PCFieldSplitRestoreSubKSP()` when the array of `KSP` is no longer needed 2520f1580f4eSBarry Smith 2521feefa0e1SJacob Faibussowitsch Developer Notes: 2522f1580f4eSBarry Smith There should be a `PCFieldSplitRestoreSubKSP()` instead of requiring the user to call `PetscFree()` 2523f1580f4eSBarry Smith 252460f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()`, `PCFieldSplitSchurGetSubKSP()` 25250971522cSBarry Smith @*/ 2526d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetSubKSP(PC pc, PetscInt *n, KSP *subksp[]) 2527d71ae5a4SJacob Faibussowitsch { 25280971522cSBarry Smith PetscFunctionBegin; 25290700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 25304f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 2531cac4c232SBarry Smith PetscUseMethod(pc, "PCFieldSplitGetSubKSP_C", (PC, PetscInt *, KSP **), (pc, n, subksp)); 25323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25330971522cSBarry Smith } 25340971522cSBarry Smith 2535285fb4e2SStefano Zampini /*@C 2536f1580f4eSBarry Smith PCFieldSplitSchurGetSubKSP - Gets the `KSP` contexts used inside the Schur complement based `PCFIELDSPLIT` 2537285fb4e2SStefano Zampini 253860f59c3bSBarry Smith Collective 2539285fb4e2SStefano Zampini 2540285fb4e2SStefano Zampini Input Parameter: 2541285fb4e2SStefano Zampini . pc - the preconditioner context 2542285fb4e2SStefano Zampini 2543285fb4e2SStefano Zampini Output Parameters: 2544285fb4e2SStefano Zampini + n - the number of splits 2545f1580f4eSBarry Smith - subksp - the array of `KSP` contexts 2546285fb4e2SStefano Zampini 2547285fb4e2SStefano Zampini Level: advanced 2548285fb4e2SStefano Zampini 2549f1580f4eSBarry Smith Notes: 2550f1580f4eSBarry Smith After `PCFieldSplitSchurGetSubKSP()` the array of `KSP`s is to be freed by the user with `PetscFree()` 2551f1580f4eSBarry Smith (not the `KSP` just the array that contains them). 2552f1580f4eSBarry Smith 2553f1580f4eSBarry Smith You must call `PCSetUp()` before calling `PCFieldSplitSchurGetSubKSP()`. 2554f1580f4eSBarry Smith 2555f1580f4eSBarry Smith If the fieldsplit type is of type `PC_COMPOSITE_SCHUR`, it returns (in order) 2556f1580f4eSBarry Smith + 1 - the `KSP` used for the (1,1) block 2557f1580f4eSBarry Smith . 2 - the `KSP` used for the Schur complement (not the one used for the interior Schur solver) 2558f1580f4eSBarry Smith - 3 - the `KSP` used for the (1,1) block in the upper triangular factor (if different from that of the (1,1) block). 2559f1580f4eSBarry Smith 2560f1580f4eSBarry Smith It returns a null array if the fieldsplit is not of type `PC_COMPOSITE_SCHUR`; in this case, you should use `PCFieldSplitGetSubKSP()`. 2561f1580f4eSBarry Smith 2562feaf08eaSBarry Smith Fortran Note: 2563e41f517fSBarry Smith Call `PCFieldSplitSchurRestoreSubKSP()` when the array of `KSP` is no longer needed 2564f1580f4eSBarry Smith 2565f1580f4eSBarry Smith Developer Notes: 2566f1580f4eSBarry Smith There should be a `PCFieldSplitRestoreSubKSP()` instead of requiring the user to call `PetscFree()` 2567f1580f4eSBarry Smith 2568f1580f4eSBarry Smith Should the functionality of `PCFieldSplitSchurGetSubKSP()` and `PCFieldSplitGetSubKSP()` be merged? 2569f1580f4eSBarry Smith 257060f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()`, `PCFieldSplitGetSubKSP()` 2571285fb4e2SStefano Zampini @*/ 2572d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSchurGetSubKSP(PC pc, PetscInt *n, KSP *subksp[]) 2573d71ae5a4SJacob Faibussowitsch { 2574285fb4e2SStefano Zampini PetscFunctionBegin; 2575285fb4e2SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 25764f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 2577cac4c232SBarry Smith PetscUseMethod(pc, "PCFieldSplitSchurGetSubKSP_C", (PC, PetscInt *, KSP **), (pc, n, subksp)); 25783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2579285fb4e2SStefano Zampini } 2580285fb4e2SStefano Zampini 2581e69d4d44SBarry Smith /*@ 2582c14e9f00SDavid Andrs PCFieldSplitSetSchurPre - Indicates from what operator the preconditioner is constructed for the Schur complement. 2583a1e3cbf9SBarry Smith The default is the A11 matrix. 2584e69d4d44SBarry Smith 2585c3339decSBarry Smith Collective 2586e69d4d44SBarry Smith 2587e69d4d44SBarry Smith Input Parameters: 2588e69d4d44SBarry Smith + pc - the preconditioner context 2589f1580f4eSBarry Smith . ptype - which matrix to use for preconditioning the Schur complement: `PC_FIELDSPLIT_SCHUR_PRE_A11` (default), 2590f1580f4eSBarry Smith `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PC_FIELDSPLIT_SCHUR_PRE_USER`, 2591f1580f4eSBarry Smith `PC_FIELDSPLIT_SCHUR_PRE_SELFP`, and `PC_FIELDSPLIT_SCHUR_PRE_FULL` 259260f59c3bSBarry Smith - pre - matrix to use for preconditioning, or `NULL` 2593084e4875SJed Brown 2594f1580f4eSBarry Smith Options Database Keys: 2595d59693daSPierre Jolivet + -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is `a11`. See notes for meaning of various arguments 2596a1e3cbf9SBarry Smith - -fieldsplit_1_pc_type <pctype> - the preconditioner algorithm that is used to construct the preconditioner from the operator 2597e69d4d44SBarry Smith 259860f59c3bSBarry Smith Level: intermediate 259960f59c3bSBarry Smith 2600fd1303e9SJungho Lee Notes: 2601f1580f4eSBarry Smith If ptype is 2602f1580f4eSBarry Smith + a11 - the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner 2603f1580f4eSBarry Smith matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix 2604f1580f4eSBarry Smith . self - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix: 2605e7593814SPierre Jolivet The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM` 2606f1580f4eSBarry Smith . user - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument 2607f1580f4eSBarry Smith to this function). 2608a077d33dSBarry Smith . selfp - the preconditioning for the Schur complement is generated from an explicitly-assembled approximation $ Sp = A11 - A10 inv(diag(A00)) A01 $ 2609f1580f4eSBarry Smith This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be 2610d59693daSPierre Jolivet lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump` 2611f1580f4eSBarry Smith - full - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation 2612f1580f4eSBarry Smith computed internally by `PCFIELDSPLIT` (this is expensive) 2613f1580f4eSBarry Smith useful mostly as a test that the Schur complement approach can work for your problem 2614fd1303e9SJungho Lee 2615d59693daSPierre Jolivet When solving a saddle point problem, where the A11 block is identically zero, using `a11` as the ptype only makes sense 2616a077d33dSBarry Smith with the additional option `-fieldsplit_1_pc_type none`. Usually for saddle point problems one would use a `ptype` of `self` and 2617d59693daSPierre Jolivet `-fieldsplit_1_pc_type lsc` which uses the least squares commutator to compute a preconditioner for the Schur complement. 2618fd1303e9SJungho Lee 2619a077d33dSBarry Smith Developer Note: 2620a077d33dSBarry Smith The name of this function and the option `-pc_fieldsplit_schur_precondition` are inconsistent; precondition should be used everywhere. 2621feefa0e1SJacob Faibussowitsch 2622a077d33dSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSchurPre()`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`, 2623a077d33dSBarry Smith `MatSchurComplementSetAinvType()`, `PCLSC`, `PCFieldSplitSetSchurFactType()` 2624e69d4d44SBarry Smith @*/ 2625d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetSchurPre(PC pc, PCFieldSplitSchurPreType ptype, Mat pre) 2626d71ae5a4SJacob Faibussowitsch { 2627e69d4d44SBarry Smith PetscFunctionBegin; 26280700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2629cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetSchurPre_C", (PC, PCFieldSplitSchurPreType, Mat), (pc, ptype, pre)); 26303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2631e69d4d44SBarry Smith } 2632686bed4dSStefano Zampini 2633d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSchurPrecondition(PC pc, PCFieldSplitSchurPreType ptype, Mat pre) 2634d71ae5a4SJacob Faibussowitsch { 26359371c9d4SSatish Balay return PCFieldSplitSetSchurPre(pc, ptype, pre); 26369371c9d4SSatish Balay } /* Deprecated name */ 2637e69d4d44SBarry Smith 263837a82bf0SJed Brown /*@ 263937a82bf0SJed Brown PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be 2640f1580f4eSBarry Smith preconditioned. See `PCFieldSplitSetSchurPre()` for details. 264137a82bf0SJed Brown 2642c3339decSBarry Smith Logically Collective 264337a82bf0SJed Brown 2644f899ff85SJose E. Roman Input Parameter: 264537a82bf0SJed Brown . pc - the preconditioner context 264637a82bf0SJed Brown 264737a82bf0SJed Brown Output Parameters: 2648d59693daSPierre Jolivet + ptype - which matrix to use for preconditioning the Schur complement: `PC_FIELDSPLIT_SCHUR_PRE_A11`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PC_FIELDSPLIT_SCHUR_PRE_USER` 2649d59693daSPierre Jolivet - pre - matrix to use for preconditioning (with `PC_FIELDSPLIT_SCHUR_PRE_USER`), or `NULL` 265037a82bf0SJed Brown 265137a82bf0SJed Brown Level: intermediate 265237a82bf0SJed Brown 2653d59693daSPierre Jolivet .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitSetSchurPre()`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`, `PCLSC` 265437a82bf0SJed Brown @*/ 2655d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetSchurPre(PC pc, PCFieldSplitSchurPreType *ptype, Mat *pre) 2656d71ae5a4SJacob Faibussowitsch { 265737a82bf0SJed Brown PetscFunctionBegin; 265837a82bf0SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2659cac4c232SBarry Smith PetscUseMethod(pc, "PCFieldSplitGetSchurPre_C", (PC, PCFieldSplitSchurPreType *, Mat *), (pc, ptype, pre)); 26603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2661e69d4d44SBarry Smith } 2662e69d4d44SBarry Smith 266345e7fc46SDmitry Karpeev /*@ 2664f1580f4eSBarry Smith PCFieldSplitSchurGetS - extract the `MATSCHURCOMPLEMENT` object used by this `PCFIELDSPLIT` in case it needs to be configured separately 266545e7fc46SDmitry Karpeev 266620f4b53cSBarry Smith Not Collective 266745e7fc46SDmitry Karpeev 266845e7fc46SDmitry Karpeev Input Parameter: 266945e7fc46SDmitry Karpeev . pc - the preconditioner context 267045e7fc46SDmitry Karpeev 267145e7fc46SDmitry Karpeev Output Parameter: 2672470b340bSDmitry Karpeev . S - the Schur complement matrix 267345e7fc46SDmitry Karpeev 267460f59c3bSBarry Smith Level: advanced 267560f59c3bSBarry Smith 2676f1580f4eSBarry Smith Note: 2677f1580f4eSBarry Smith This matrix should not be destroyed using `MatDestroy()`; rather, use `PCFieldSplitSchurRestoreS()`. 267845e7fc46SDmitry Karpeev 267960f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSchurPreType`, `PCFieldSplitSetSchurPre()`, `MATSCHURCOMPLEMENT`, `PCFieldSplitSchurRestoreS()`, 2680f1580f4eSBarry Smith `MatCreateSchurComplement()`, `MatSchurComplementGetKSP()`, `MatSchurComplementComputeExplicitOperator()`, `MatGetSchurComplement()` 268145e7fc46SDmitry Karpeev @*/ 2682d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSchurGetS(PC pc, Mat *S) 2683d71ae5a4SJacob Faibussowitsch { 268445e7fc46SDmitry Karpeev const char *t; 268545e7fc46SDmitry Karpeev PetscBool isfs; 268645e7fc46SDmitry Karpeev PC_FieldSplit *jac; 268745e7fc46SDmitry Karpeev 268845e7fc46SDmitry Karpeev PetscFunctionBegin; 268945e7fc46SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 26909566063dSJacob Faibussowitsch PetscCall(PetscObjectGetType((PetscObject)pc, &t)); 26919566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(t, PCFIELDSPLIT, &isfs)); 269228b400f6SJacob Faibussowitsch PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PC of type PCFIELDSPLIT, got %s instead", t); 269345e7fc46SDmitry Karpeev jac = (PC_FieldSplit *)pc->data; 269463a3b9bcSJacob Faibussowitsch PetscCheck(jac->type == PC_COMPOSITE_SCHUR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PCFIELDSPLIT of type SCHUR, got %d instead", jac->type); 2695470b340bSDmitry Karpeev if (S) *S = jac->schur; 26963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269745e7fc46SDmitry Karpeev } 269845e7fc46SDmitry Karpeev 2699470b340bSDmitry Karpeev /*@ 2700f1580f4eSBarry Smith PCFieldSplitSchurRestoreS - returns the `MATSCHURCOMPLEMENT` matrix used by this `PC` 2701470b340bSDmitry Karpeev 270220f4b53cSBarry Smith Not Collective 2703470b340bSDmitry Karpeev 2704470b340bSDmitry Karpeev Input Parameters: 2705470b340bSDmitry Karpeev + pc - the preconditioner context 2706a2b725a8SWilliam Gropp - S - the Schur complement matrix 2707470b340bSDmitry Karpeev 2708470b340bSDmitry Karpeev Level: advanced 2709470b340bSDmitry Karpeev 271060f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSchurPreType`, `PCFieldSplitSetSchurPre()`, `MatSchurComplement`, `PCFieldSplitSchurGetS()` 2711470b340bSDmitry Karpeev @*/ 2712d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSchurRestoreS(PC pc, Mat *S) 2713d71ae5a4SJacob Faibussowitsch { 2714470b340bSDmitry Karpeev const char *t; 2715470b340bSDmitry Karpeev PetscBool isfs; 2716470b340bSDmitry Karpeev PC_FieldSplit *jac; 2717470b340bSDmitry Karpeev 2718470b340bSDmitry Karpeev PetscFunctionBegin; 2719470b340bSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 27209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetType((PetscObject)pc, &t)); 27219566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(t, PCFIELDSPLIT, &isfs)); 272228b400f6SJacob Faibussowitsch PetscCheck(isfs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PC of type PCFIELDSPLIT, got %s instead", t); 2723470b340bSDmitry Karpeev jac = (PC_FieldSplit *)pc->data; 272463a3b9bcSJacob Faibussowitsch PetscCheck(jac->type == PC_COMPOSITE_SCHUR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected PCFIELDSPLIT of type SCHUR, got %d instead", jac->type); 27257827d75bSBarry Smith PetscCheck(S && (*S == jac->schur), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatSchurComplement restored is not the same as gotten"); 27263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2727470b340bSDmitry Karpeev } 2728470b340bSDmitry Karpeev 2729d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetSchurPre_FieldSplit(PC pc, PCFieldSplitSchurPreType ptype, Mat pre) 2730d71ae5a4SJacob Faibussowitsch { 2731e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2732e69d4d44SBarry Smith 2733e69d4d44SBarry Smith PetscFunctionBegin; 2734084e4875SJed Brown jac->schurpre = ptype; 2735a7476a74SDmitry Karpeev if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) { 27369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->schur_user)); 2737084e4875SJed Brown jac->schur_user = pre; 27389566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)jac->schur_user)); 2739084e4875SJed Brown } 27403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2741e69d4d44SBarry Smith } 2742e69d4d44SBarry Smith 2743d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitGetSchurPre_FieldSplit(PC pc, PCFieldSplitSchurPreType *ptype, Mat *pre) 2744d71ae5a4SJacob Faibussowitsch { 274537a82bf0SJed Brown PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 274637a82bf0SJed Brown 274737a82bf0SJed Brown PetscFunctionBegin; 27486056e507SPierre Jolivet if (ptype) *ptype = jac->schurpre; 27496056e507SPierre Jolivet if (pre) *pre = jac->schur_user; 27503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275137a82bf0SJed Brown } 275237a82bf0SJed Brown 2753ab1df9f5SJed Brown /*@ 27541d27aa22SBarry Smith PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain in the preconditioner {cite}`murphy2000note` and {cite}`ipsen2001note` 2755ab1df9f5SJed Brown 2756c3339decSBarry Smith Collective 2757ab1df9f5SJed Brown 2758ab1df9f5SJed Brown Input Parameters: 2759ab1df9f5SJed Brown + pc - the preconditioner context 2760f1580f4eSBarry Smith - ftype - which blocks of factorization to retain, `PC_FIELDSPLIT_SCHUR_FACT_FULL` is default 2761ab1df9f5SJed Brown 2762f1580f4eSBarry Smith Options Database Key: 2763d59693daSPierre Jolivet . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - default is `full` 2764ab1df9f5SJed Brown 2765ab1df9f5SJed Brown Level: intermediate 2766ab1df9f5SJed Brown 2767ab1df9f5SJed Brown Notes: 2768366f9a6aSPierre Jolivet The `full` factorization is 2769ab1df9f5SJed Brown 2770a077d33dSBarry Smith ```{math} 2771a077d33dSBarry Smith \left(\begin{array}{cc} A & B \\ 2772a077d33dSBarry Smith C & E \\ 2773a077d33dSBarry Smith \end{array}\right) = 2774366f9a6aSPierre Jolivet \left(\begin{array}{cc} I & 0 \\ 2775366f9a6aSPierre Jolivet C A^{-1} & I \\ 2776a077d33dSBarry Smith \end{array}\right) 2777a077d33dSBarry Smith \left(\begin{array}{cc} A & 0 \\ 2778a077d33dSBarry Smith 0 & S \\ 2779a077d33dSBarry Smith \end{array}\right) 2780a077d33dSBarry Smith \left(\begin{array}{cc} I & A^{-1}B \\ 2781a077d33dSBarry Smith 0 & I \\ 2782366f9a6aSPierre Jolivet \end{array}\right) = L D U, 2783a077d33dSBarry Smith ``` 2784a077d33dSBarry Smith 2785366f9a6aSPierre Jolivet where $ S = E - C A^{-1} B $. In practice, the full factorization is applied via block triangular solves with the grouping $L(DU)$. `upper` uses $DU$, `lower` uses $LD$, 2786366f9a6aSPierre Jolivet and `diag` is the diagonal part with the sign of $S$ flipped (because this makes the preconditioner positive definite for many formulations, 2787a077d33dSBarry Smith thus allowing the use of `KSPMINRES)`. Sign flipping of $S$ can be turned off with `PCFieldSplitSetSchurScale()`. 2788a077d33dSBarry Smith 2789a077d33dSBarry Smith If $A$ and $S$ are solved exactly 2790366f9a6aSPierre Jolivet + 1 - `full` factorization is a direct solver. 2791366f9a6aSPierre Jolivet . 2 - 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. 2792366f9a6aSPierre Jolivet - 3 - 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. 2793ab1df9f5SJed Brown 2794f1580f4eSBarry Smith If the iteration count is very low, consider using `KSPFGMRES` or `KSPGCR` which can use one less preconditioner 27950ffb0e17SBarry 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. 27960ffb0e17SBarry Smith 2797366f9a6aSPierre Jolivet For symmetric problems in which $A$ is positive definite and $S$ is negative definite, `diag` can be used with `KSPMINRES`. 27980ffb0e17SBarry Smith 2799366f9a6aSPierre Jolivet A flexible method like `KSPFGMRES` or `KSPGCR`, [](sec_flexibleksp), 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$). 2800ab1df9f5SJed Brown 28011d27aa22SBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFieldSplitGetSubKSP()`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurPreType`, `PCFieldSplitSetSchurScale()`, 2802a077d33dSBarry Smith [](sec_flexibleksp), `PCFieldSplitSetSchurPre()` 2803ab1df9f5SJed Brown @*/ 2804d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetSchurFactType(PC pc, PCFieldSplitSchurFactType ftype) 2805d71ae5a4SJacob Faibussowitsch { 2806ab1df9f5SJed Brown PetscFunctionBegin; 2807ab1df9f5SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2808cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetSchurFactType_C", (PC, PCFieldSplitSchurFactType), (pc, ftype)); 28093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2810ab1df9f5SJed Brown } 2811ab1df9f5SJed Brown 2812d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc, PCFieldSplitSchurFactType ftype) 2813d71ae5a4SJacob Faibussowitsch { 2814ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2815ab1df9f5SJed Brown 2816ab1df9f5SJed Brown PetscFunctionBegin; 2817ab1df9f5SJed Brown jac->schurfactorization = ftype; 28183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2819ab1df9f5SJed Brown } 2820ab1df9f5SJed Brown 2821c096484dSStefano Zampini /*@ 2822f1580f4eSBarry Smith PCFieldSplitSetSchurScale - Controls the sign flip of S for `PC_FIELDSPLIT_SCHUR_FACT_DIAG`. 2823c096484dSStefano Zampini 2824c3339decSBarry Smith Collective 2825c096484dSStefano Zampini 2826c096484dSStefano Zampini Input Parameters: 2827c096484dSStefano Zampini + pc - the preconditioner context 2828c096484dSStefano Zampini - scale - scaling factor for the Schur complement 2829c096484dSStefano Zampini 2830f1580f4eSBarry Smith Options Database Key: 28311d27aa22SBarry Smith . -pc_fieldsplit_schur_scale <scale> - default is -1.0 2832c096484dSStefano Zampini 2833c096484dSStefano Zampini Level: intermediate 2834c096484dSStefano Zampini 283542747ad1SJacob Faibussowitsch .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetFields()`, `PCFieldSplitSchurFactType`, `PCFieldSplitSetSchurFactType()` 2836c096484dSStefano Zampini @*/ 2837d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetSchurScale(PC pc, PetscScalar scale) 2838d71ae5a4SJacob Faibussowitsch { 2839c096484dSStefano Zampini PetscFunctionBegin; 2840c096484dSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2841c096484dSStefano Zampini PetscValidLogicalCollectiveScalar(pc, scale, 2); 2842cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetSchurScale_C", (PC, PetscScalar), (pc, scale)); 28433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2844c096484dSStefano Zampini } 2845c096484dSStefano Zampini 2846d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc, PetscScalar scale) 2847d71ae5a4SJacob Faibussowitsch { 2848c096484dSStefano Zampini PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2849c096484dSStefano Zampini 2850c096484dSStefano Zampini PetscFunctionBegin; 2851c096484dSStefano Zampini jac->schurscale = scale; 28523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2853c096484dSStefano Zampini } 2854c096484dSStefano Zampini 285530ad9308SMatthew Knepley /*@C 28568c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 285730ad9308SMatthew Knepley 2858c3339decSBarry Smith Collective 285930ad9308SMatthew Knepley 286030ad9308SMatthew Knepley Input Parameter: 286130ad9308SMatthew Knepley . pc - the preconditioner context 286230ad9308SMatthew Knepley 286330ad9308SMatthew Knepley Output Parameters: 2864a04f6461SBarry Smith + A00 - the (0,0) block 2865a04f6461SBarry Smith . A01 - the (0,1) block 2866a04f6461SBarry Smith . A10 - the (1,0) block 2867a04f6461SBarry Smith - A11 - the (1,1) block 286830ad9308SMatthew Knepley 286930ad9308SMatthew Knepley Level: advanced 287030ad9308SMatthew Knepley 28715d83a8b1SBarry Smith Note: 28725d83a8b1SBarry Smith Use `NULL` for any unneeded output arguments 28735d83a8b1SBarry Smith 287460f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `MatSchurComplementGetSubMatrices()`, `MatSchurComplementSetSubMatrices()` 287530ad9308SMatthew Knepley @*/ 2876d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc, Mat *A00, Mat *A01, Mat *A10, Mat *A11) 2877d71ae5a4SJacob Faibussowitsch { 287830ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 287930ad9308SMatthew Knepley 288030ad9308SMatthew Knepley PetscFunctionBegin; 28810700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 288208401ef6SPierre Jolivet PetscCheck(jac->type == PC_COMPOSITE_SCHUR, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 2883a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 2884a04f6461SBarry Smith if (A01) *A01 = jac->B; 2885a04f6461SBarry Smith if (A10) *A10 = jac->C; 2886a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 28873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 288830ad9308SMatthew Knepley } 288930ad9308SMatthew Knepley 2890b09e027eSCarola Kruse /*@ 28911d27aa22SBarry Smith PCFieldSplitSetGKBTol - Sets the solver tolerance for the generalized Golub-Kahan bidiagonalization preconditioner {cite}`arioli2013` in `PCFIELDSPLIT` 2892b09e027eSCarola Kruse 2893c3339decSBarry Smith Collective 2894e071a0a4SCarola Kruse 2895b09e027eSCarola Kruse Input Parameters: 2896b09e027eSCarola Kruse + pc - the preconditioner context 2897b09e027eSCarola Kruse - tolerance - the solver tolerance 2898b09e027eSCarola Kruse 2899f1580f4eSBarry Smith Options Database Key: 29001d27aa22SBarry Smith . -pc_fieldsplit_gkb_tol <tolerance> - default is 1e-5 2901b09e027eSCarola Kruse 2902b09e027eSCarola Kruse Level: intermediate 2903b09e027eSCarola Kruse 2904f1580f4eSBarry Smith Note: 29051d27aa22SBarry Smith The generalized GKB algorithm {cite}`arioli2013` uses a lower bound estimate of the error in energy norm as stopping criterion. 2906f1580f4eSBarry Smith It stops once the lower bound estimate undershoots the required solver tolerance. Although the actual error might be bigger than 29071d27aa22SBarry Smith this estimate, the stopping criterion is satisfactory in practical cases. 2908f1580f4eSBarry Smith 290960f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetGKBDelay()`, `PCFieldSplitSetGKBNu()`, `PCFieldSplitSetGKBMaxit()` 2910b09e027eSCarola Kruse @*/ 2911d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetGKBTol(PC pc, PetscReal tolerance) 2912d71ae5a4SJacob Faibussowitsch { 2913b09e027eSCarola Kruse PetscFunctionBegin; 2914b09e027eSCarola Kruse PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2915de482cd7SCarola Kruse PetscValidLogicalCollectiveReal(pc, tolerance, 2); 2916cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetGKBTol_C", (PC, PetscReal), (pc, tolerance)); 29173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2918b09e027eSCarola Kruse } 2919b09e027eSCarola Kruse 2920d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetGKBTol_FieldSplit(PC pc, PetscReal tolerance) 2921d71ae5a4SJacob Faibussowitsch { 2922b09e027eSCarola Kruse PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2923b09e027eSCarola Kruse 2924b09e027eSCarola Kruse PetscFunctionBegin; 2925b09e027eSCarola Kruse jac->gkbtol = tolerance; 29263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2927b09e027eSCarola Kruse } 2928b09e027eSCarola Kruse 2929b09e027eSCarola Kruse /*@ 293080670ca5SBarry Smith PCFieldSplitSetGKBMaxit - Sets the maximum number of iterations for the generalized Golub-Kahan bidiagonalization preconditioner {cite}`arioli2013` in `PCFIELDSPLIT` 2931b09e027eSCarola Kruse 2932c3339decSBarry Smith Collective 2933b09e027eSCarola Kruse 2934b09e027eSCarola Kruse Input Parameters: 2935b09e027eSCarola Kruse + pc - the preconditioner context 2936b09e027eSCarola Kruse - maxit - the maximum number of iterations 2937b09e027eSCarola Kruse 2938f1580f4eSBarry Smith Options Database Key: 29391d27aa22SBarry Smith . -pc_fieldsplit_gkb_maxit <maxit> - default is 100 2940b09e027eSCarola Kruse 2941b09e027eSCarola Kruse Level: intermediate 2942b09e027eSCarola Kruse 294360f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetGKBDelay()`, `PCFieldSplitSetGKBTol()`, `PCFieldSplitSetGKBNu()` 2944b09e027eSCarola Kruse @*/ 2945d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetGKBMaxit(PC pc, PetscInt maxit) 2946d71ae5a4SJacob Faibussowitsch { 2947b09e027eSCarola Kruse PetscFunctionBegin; 2948b09e027eSCarola Kruse PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2949de482cd7SCarola Kruse PetscValidLogicalCollectiveInt(pc, maxit, 2); 2950cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetGKBMaxit_C", (PC, PetscInt), (pc, maxit)); 29513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2952b09e027eSCarola Kruse } 2953b09e027eSCarola Kruse 2954d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetGKBMaxit_FieldSplit(PC pc, PetscInt maxit) 2955d71ae5a4SJacob Faibussowitsch { 2956b09e027eSCarola Kruse PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2957b09e027eSCarola Kruse 2958b09e027eSCarola Kruse PetscFunctionBegin; 2959b09e027eSCarola Kruse jac->gkbmaxit = maxit; 29603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2961b09e027eSCarola Kruse } 2962b09e027eSCarola Kruse 2963b09e027eSCarola Kruse /*@ 29641d27aa22SBarry Smith PCFieldSplitSetGKBDelay - Sets the delay in the lower bound error estimate in the generalized Golub-Kahan bidiagonalization {cite}`arioli2013` in `PCFIELDSPLIT` 2965e071a0a4SCarola Kruse preconditioner. 2966b09e027eSCarola Kruse 2967c3339decSBarry Smith Collective 2968b09e027eSCarola Kruse 2969b09e027eSCarola Kruse Input Parameters: 2970b09e027eSCarola Kruse + pc - the preconditioner context 2971b09e027eSCarola Kruse - delay - the delay window in the lower bound estimate 2972b09e027eSCarola Kruse 2973f1580f4eSBarry Smith Options Database Key: 29741d27aa22SBarry Smith . -pc_fieldsplit_gkb_delay <delay> - default is 5 2975b09e027eSCarola Kruse 2976b09e027eSCarola Kruse Level: intermediate 2977b09e027eSCarola Kruse 29781d27aa22SBarry Smith Notes: 29791d27aa22SBarry 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 $ 29801d27aa22SBarry 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 29811d27aa22SBarry Smith at least (`delay` + 1) iterations to stop. 2982f1580f4eSBarry Smith 29831d27aa22SBarry Smith For more details on the generalized Golub-Kahan bidiagonalization method and its lower bound stopping criterion, please refer to {cite}`arioli2013` 2984f1580f4eSBarry Smith 298560f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetGKBNu()`, `PCFieldSplitSetGKBTol()`, `PCFieldSplitSetGKBMaxit()` 2986b09e027eSCarola Kruse @*/ 2987d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetGKBDelay(PC pc, PetscInt delay) 2988d71ae5a4SJacob Faibussowitsch { 2989b09e027eSCarola Kruse PetscFunctionBegin; 2990b09e027eSCarola Kruse PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2991b09e027eSCarola Kruse PetscValidLogicalCollectiveInt(pc, delay, 2); 2992cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetGKBDelay_C", (PC, PetscInt), (pc, delay)); 29933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2994b09e027eSCarola Kruse } 2995b09e027eSCarola Kruse 2996d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetGKBDelay_FieldSplit(PC pc, PetscInt delay) 2997d71ae5a4SJacob Faibussowitsch { 2998b09e027eSCarola Kruse PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 2999b09e027eSCarola Kruse 3000b09e027eSCarola Kruse PetscFunctionBegin; 3001b09e027eSCarola Kruse jac->gkbdelay = delay; 30023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3003b09e027eSCarola Kruse } 3004b09e027eSCarola Kruse 3005b09e027eSCarola Kruse /*@ 30061d27aa22SBarry Smith PCFieldSplitSetGKBNu - Sets the scalar value nu >= 0 in the transformation H = A00 + nu*A01*A01' of the (1,1) block in the 30071d27aa22SBarry Smith Golub-Kahan bidiagonalization preconditioner {cite}`arioli2013` in `PCFIELDSPLIT` 3008b09e027eSCarola Kruse 3009c3339decSBarry Smith Collective 3010f1580f4eSBarry Smith 3011f1580f4eSBarry Smith Input Parameters: 3012f1580f4eSBarry Smith + pc - the preconditioner context 3013f1580f4eSBarry Smith - nu - the shift parameter 3014f1580f4eSBarry Smith 301520f4b53cSBarry Smith Options Database Key: 30161d27aa22SBarry Smith . -pc_fieldsplit_gkb_nu <nu> - default is 1 3017f1580f4eSBarry Smith 3018f1580f4eSBarry Smith Level: intermediate 3019b09e027eSCarola Kruse 3020b09e027eSCarola Kruse Notes: 30211d27aa22SBarry Smith This shift is in general done to obtain better convergence properties for the outer loop of the algorithm. This is often achieved by choosing `nu` sufficiently large. However, 30221d27aa22SBarry Smith if `nu` is chosen too large, the matrix H might be badly conditioned and the solution of the linear system $Hx = b$ in the inner loop becomes difficult. It is therefore 3023b09e027eSCarola Kruse necessary to find a good balance in between the convergence of the inner and outer loop. 3024b09e027eSCarola Kruse 30251d27aa22SBarry Smith For `nu` = 0, no shift is done. In this case A00 has to be positive definite. The matrix N in {cite}`arioli2013` is then chosen as identity. 3026b09e027eSCarola Kruse 302760f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetGKBDelay()`, `PCFieldSplitSetGKBTol()`, `PCFieldSplitSetGKBMaxit()` 3028b09e027eSCarola Kruse @*/ 3029d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetGKBNu(PC pc, PetscReal nu) 3030d71ae5a4SJacob Faibussowitsch { 3031b09e027eSCarola Kruse PetscFunctionBegin; 3032b09e027eSCarola Kruse PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3033de482cd7SCarola Kruse PetscValidLogicalCollectiveReal(pc, nu, 2); 3034cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetGKBNu_C", (PC, PetscReal), (pc, nu)); 30353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3036b09e027eSCarola Kruse } 3037b09e027eSCarola Kruse 3038d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetGKBNu_FieldSplit(PC pc, PetscReal nu) 3039d71ae5a4SJacob Faibussowitsch { 3040b09e027eSCarola Kruse PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 3041b09e027eSCarola Kruse 3042b09e027eSCarola Kruse PetscFunctionBegin; 3043b09e027eSCarola Kruse jac->gkbnu = nu; 30443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3045b09e027eSCarola Kruse } 3046b09e027eSCarola Kruse 3047d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc, PCCompositeType type) 3048d71ae5a4SJacob Faibussowitsch { 304979416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 305079416396SBarry Smith 305179416396SBarry Smith PetscFunctionBegin; 305279416396SBarry Smith jac->type = type; 30532e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", NULL)); 30542e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurPre_C", NULL)); 30552e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSchurPre_C", NULL)); 30562e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurFactType_C", NULL)); 30572e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurScale_C", NULL)); 30582e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBTol_C", NULL)); 30592e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBMaxit_C", NULL)); 30602e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBNu_C", NULL)); 30612e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBDelay_C", NULL)); 3062e071a0a4SCarola Kruse 30633b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 30643b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 30657b665727SPierre Jolivet pc->ops->applytranspose = PCApplyTranspose_FieldSplit_Schur; 30663b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 306773716367SStefano Zampini pc->ops->setuponblocks = PCSetUpOnBlocks_FieldSplit_Schur; 30682fa5cd67SKarl Rupp 30699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit_Schur)); 30709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurPre_C", PCFieldSplitSetSchurPre_FieldSplit)); 30719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSchurPre_C", PCFieldSplitGetSchurPre_FieldSplit)); 30729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurFactType_C", PCFieldSplitSetSchurFactType_FieldSplit)); 30739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetSchurScale_C", PCFieldSplitSetSchurScale_FieldSplit)); 3074a51937d4SCarola Kruse } else if (type == PC_COMPOSITE_GKB) { 3075a51937d4SCarola Kruse pc->ops->apply = PCApply_FieldSplit_GKB; 30767ff38633SStefano Zampini pc->ops->applytranspose = NULL; 3077a51937d4SCarola Kruse pc->ops->view = PCView_FieldSplit_GKB; 307873716367SStefano Zampini pc->ops->setuponblocks = PCSetUpOnBlocks_FieldSplit_GKB; 3079e69d4d44SBarry Smith 30809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit)); 30819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBTol_C", PCFieldSplitSetGKBTol_FieldSplit)); 30829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBMaxit_C", PCFieldSplitSetGKBMaxit_FieldSplit)); 30839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBNu_C", PCFieldSplitSetGKBNu_FieldSplit)); 30849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetGKBDelay_C", PCFieldSplitSetGKBDelay_FieldSplit)); 30853b224e63SBarry Smith } else { 30863b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 30877ff38633SStefano Zampini pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 30883b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 308973716367SStefano Zampini pc->ops->setuponblocks = PCSetUpOnBlocks_FieldSplit; 30902fa5cd67SKarl Rupp 30919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitGetSubKSP_C", PCFieldSplitGetSubKSP_FieldSplit)); 30923b224e63SBarry Smith } 30933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 309479416396SBarry Smith } 309579416396SBarry Smith 3096d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc, PetscInt bs) 3097d71ae5a4SJacob Faibussowitsch { 309851f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 309951f519a2SBarry Smith 310051f519a2SBarry Smith PetscFunctionBegin; 310163a3b9bcSJacob Faibussowitsch PetscCheck(bs >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %" PetscInt_FMT, bs); 31022472a847SBarry 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); 310351f519a2SBarry Smith jac->bs = bs; 31043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310551f519a2SBarry Smith } 310651f519a2SBarry Smith 3107d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_FieldSplit(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[]) 3108d71ae5a4SJacob Faibussowitsch { 31095ddf11f8SNicolas Barnafi PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 31105ddf11f8SNicolas Barnafi PC_FieldSplitLink ilink_current = jac->head; 31115ddf11f8SNicolas Barnafi IS is_owned; 31125ddf11f8SNicolas Barnafi 31135ddf11f8SNicolas Barnafi PetscFunctionBegin; 31145ddf11f8SNicolas Barnafi jac->coordinates_set = PETSC_TRUE; // Internal flag 3115f3fa974cSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(pc->mat, &is_owned, NULL)); 31165ddf11f8SNicolas Barnafi 31175ddf11f8SNicolas Barnafi while (ilink_current) { 31185ddf11f8SNicolas Barnafi // For each IS, embed it to get local coords indces 31195ddf11f8SNicolas Barnafi IS is_coords; 31205ddf11f8SNicolas Barnafi PetscInt ndofs_block; 31215ddf11f8SNicolas Barnafi const PetscInt *block_dofs_enumeration; // Numbering of the dofs relevant to the current block 31225ddf11f8SNicolas Barnafi 31235ddf11f8SNicolas Barnafi // Setting drop to true for safety. It should make no difference. 31249566063dSJacob Faibussowitsch PetscCall(ISEmbed(ilink_current->is, is_owned, PETSC_TRUE, &is_coords)); 31259566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is_coords, &ndofs_block)); 31269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_coords, &block_dofs_enumeration)); 31275ddf11f8SNicolas Barnafi 31285ddf11f8SNicolas Barnafi // Allocate coordinates vector and set it directly 3129f4f49eeaSPierre Jolivet PetscCall(PetscMalloc1(ndofs_block * dim, &ilink_current->coords)); 31305ddf11f8SNicolas Barnafi for (PetscInt dof = 0; dof < ndofs_block; ++dof) { 3131ad540459SPierre Jolivet for (PetscInt d = 0; d < dim; ++d) (ilink_current->coords)[dim * dof + d] = coords[dim * block_dofs_enumeration[dof] + d]; 31325ddf11f8SNicolas Barnafi } 31335ddf11f8SNicolas Barnafi ilink_current->dim = dim; 31345ddf11f8SNicolas Barnafi ilink_current->ndofs = ndofs_block; 31359566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_coords, &block_dofs_enumeration)); 31369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_coords)); 31375ddf11f8SNicolas Barnafi ilink_current = ilink_current->next; 31385ddf11f8SNicolas Barnafi } 31399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_owned)); 31403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31415ddf11f8SNicolas Barnafi } 31425ddf11f8SNicolas Barnafi 3143bc08b0f1SBarry Smith /*@ 3144f1580f4eSBarry Smith PCFieldSplitSetType - Sets the type, `PCCompositeType`, of a `PCFIELDSPLIT` 314579416396SBarry Smith 3146c3339decSBarry Smith Collective 314779416396SBarry Smith 3148d8d19677SJose E. Roman Input Parameters: 3149a2b725a8SWilliam Gropp + pc - the preconditioner context 3150a077d33dSBarry Smith - type - `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE` (default), `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR`, 3151a077d33dSBarry Smith `PC_COMPOSITE_GKB` 315279416396SBarry Smith 315379416396SBarry Smith Options Database Key: 31541d27aa22SBarry Smith . -pc_fieldsplit_type <one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 315579416396SBarry Smith 3156feefa0e1SJacob Faibussowitsch Level: intermediate 315779416396SBarry Smith 315860f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCCompositeType`, `PCCompositeGetType()`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, 3159a077d33dSBarry Smith `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR`, `PCFieldSplitSetSchurFactType()` 316079416396SBarry Smith @*/ 3161d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetType(PC pc, PCCompositeType type) 3162d71ae5a4SJacob Faibussowitsch { 316379416396SBarry Smith PetscFunctionBegin; 31640700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3165cac4c232SBarry Smith PetscTryMethod(pc, "PCFieldSplitSetType_C", (PC, PCCompositeType), (pc, type)); 31663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316779416396SBarry Smith } 316879416396SBarry Smith 3169b02e2d75SMatthew G Knepley /*@ 3170f1580f4eSBarry Smith PCFieldSplitGetType - Gets the type, `PCCompositeType`, of a `PCFIELDSPLIT` 3171b02e2d75SMatthew G Knepley 3172b02e2d75SMatthew G Knepley Not collective 3173b02e2d75SMatthew G Knepley 3174b02e2d75SMatthew G Knepley Input Parameter: 3175b02e2d75SMatthew G Knepley . pc - the preconditioner context 3176b02e2d75SMatthew G Knepley 3177b02e2d75SMatthew G Knepley Output Parameter: 3178f1580f4eSBarry Smith . type - `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE` (default), `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR` 3179b02e2d75SMatthew G Knepley 3180feefa0e1SJacob Faibussowitsch Level: intermediate 3181b02e2d75SMatthew G Knepley 318260f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCCompositeSetType()`, `PCFIELDSPLIT`, `PCCompositeType`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, 3183f1580f4eSBarry Smith `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PC_COMPOSITE_SCHUR` 3184b02e2d75SMatthew G Knepley @*/ 3185d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 3186d71ae5a4SJacob Faibussowitsch { 3187b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 3188b02e2d75SMatthew G Knepley 3189b02e2d75SMatthew G Knepley PetscFunctionBegin; 3190b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 31914f572ea9SToby Isaac PetscAssertPointer(type, 2); 3192b02e2d75SMatthew G Knepley *type = jac->type; 31933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3194b02e2d75SMatthew G Knepley } 3195b02e2d75SMatthew G Knepley 31964ab8060aSDmitry Karpeev /*@ 3197f1580f4eSBarry Smith PCFieldSplitSetDMSplits - Flags whether `DMCreateFieldDecomposition()` should be used to define the splits in a `PCFIELDSPLIT`, whenever possible. 31984ab8060aSDmitry Karpeev 3199c3339decSBarry Smith Logically Collective 32004ab8060aSDmitry Karpeev 32014ab8060aSDmitry Karpeev Input Parameters: 32024ab8060aSDmitry Karpeev + pc - the preconditioner context 3203f1580f4eSBarry Smith - flg - boolean indicating whether to use field splits defined by the `DM` 32044ab8060aSDmitry Karpeev 32054ab8060aSDmitry Karpeev Options Database Key: 3206f1580f4eSBarry Smith . -pc_fieldsplit_dm_splits <bool> - use the field splits defined by the `DM` 32074ab8060aSDmitry Karpeev 3208feefa0e1SJacob Faibussowitsch Level: intermediate 32094ab8060aSDmitry Karpeev 321073ff1848SBarry Smith Developer Note: 321173ff1848SBarry Smith The name should be `PCFieldSplitSetUseDMSplits()`, similar change to options database 321273ff1848SBarry Smith 3213f8d70eaaSPierre Jolivet .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitGetDMSplits()`, `DMCreateFieldDecomposition()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()` 32144ab8060aSDmitry Karpeev @*/ 3215d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetDMSplits(PC pc, PetscBool flg) 3216d71ae5a4SJacob Faibussowitsch { 32174ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 32184ab8060aSDmitry Karpeev PetscBool isfs; 32194ab8060aSDmitry Karpeev 32204ab8060aSDmitry Karpeev PetscFunctionBegin; 32214ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 32224ab8060aSDmitry Karpeev PetscValidLogicalCollectiveBool(pc, flg, 2); 32239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs)); 3224ad540459SPierre Jolivet if (isfs) jac->dm_splits = flg; 32253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32264ab8060aSDmitry Karpeev } 32274ab8060aSDmitry Karpeev 32284ab8060aSDmitry Karpeev /*@ 3229f1580f4eSBarry Smith PCFieldSplitGetDMSplits - Returns flag indicating whether `DMCreateFieldDecomposition()` should be used to define the splits in a `PCFIELDSPLIT`, whenever possible. 32304ab8060aSDmitry Karpeev 32314ab8060aSDmitry Karpeev Logically Collective 32324ab8060aSDmitry Karpeev 32334ab8060aSDmitry Karpeev Input Parameter: 32344ab8060aSDmitry Karpeev . pc - the preconditioner context 32354ab8060aSDmitry Karpeev 32364ab8060aSDmitry Karpeev Output Parameter: 3237f1580f4eSBarry Smith . flg - boolean indicating whether to use field splits defined by the `DM` 32384ab8060aSDmitry Karpeev 3239feefa0e1SJacob Faibussowitsch Level: intermediate 32404ab8060aSDmitry Karpeev 324173ff1848SBarry Smith Developer Note: 324273ff1848SBarry Smith The name should be `PCFieldSplitGetUseDMSplits()` 324373ff1848SBarry Smith 3244f8d70eaaSPierre Jolivet .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetDMSplits()`, `DMCreateFieldDecomposition()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetIS()` 32454ab8060aSDmitry Karpeev @*/ 3246d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetDMSplits(PC pc, PetscBool *flg) 3247d71ae5a4SJacob Faibussowitsch { 32484ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 32494ab8060aSDmitry Karpeev PetscBool isfs; 32504ab8060aSDmitry Karpeev 32514ab8060aSDmitry Karpeev PetscFunctionBegin; 32524ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 32534f572ea9SToby Isaac PetscAssertPointer(flg, 2); 32549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &isfs)); 32554ab8060aSDmitry Karpeev if (isfs) { 32564ab8060aSDmitry Karpeev if (flg) *flg = jac->dm_splits; 32574ab8060aSDmitry Karpeev } 32583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32594ab8060aSDmitry Karpeev } 32604ab8060aSDmitry Karpeev 32617b752e3dSPatrick Sanan /*@ 3262f1580f4eSBarry Smith PCFieldSplitGetDetectSaddlePoint - Returns flag indicating whether `PCFIELDSPLIT` will attempt to automatically determine fields based on zero diagonal entries. 32637b752e3dSPatrick Sanan 32647b752e3dSPatrick Sanan Logically Collective 32657b752e3dSPatrick Sanan 32667b752e3dSPatrick Sanan Input Parameter: 32677b752e3dSPatrick Sanan . pc - the preconditioner context 32687b752e3dSPatrick Sanan 32697b752e3dSPatrick Sanan Output Parameter: 32707b752e3dSPatrick Sanan . flg - boolean indicating whether to detect fields or not 32717b752e3dSPatrick Sanan 3272feefa0e1SJacob Faibussowitsch Level: intermediate 32737b752e3dSPatrick Sanan 327460f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitSetDetectSaddlePoint()` 32757b752e3dSPatrick Sanan @*/ 3276d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC pc, PetscBool *flg) 3277d71ae5a4SJacob Faibussowitsch { 32787b752e3dSPatrick Sanan PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 32797b752e3dSPatrick Sanan 32807b752e3dSPatrick Sanan PetscFunctionBegin; 32817b752e3dSPatrick Sanan *flg = jac->detect; 32823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32837b752e3dSPatrick Sanan } 32847b752e3dSPatrick Sanan 32857b752e3dSPatrick Sanan /*@ 3286f1580f4eSBarry Smith PCFieldSplitSetDetectSaddlePoint - Sets flag indicating whether `PCFIELDSPLIT` will attempt to automatically determine fields based on zero diagonal entries. 32877b752e3dSPatrick Sanan 32887b752e3dSPatrick Sanan Logically Collective 32897b752e3dSPatrick Sanan 32907b752e3dSPatrick Sanan Input Parameter: 32917b752e3dSPatrick Sanan . pc - the preconditioner context 32927b752e3dSPatrick Sanan 32937b752e3dSPatrick Sanan Output Parameter: 32947b752e3dSPatrick Sanan . flg - boolean indicating whether to detect fields or not 32957b752e3dSPatrick Sanan 32967b752e3dSPatrick Sanan Options Database Key: 3297147403d9SBarry Smith . -pc_fieldsplit_detect_saddle_point <bool> - detect and use the saddle point 3298147403d9SBarry Smith 3299feefa0e1SJacob Faibussowitsch Level: intermediate 330060f59c3bSBarry Smith 3301f1580f4eSBarry Smith Note: 3302f1580f4eSBarry 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()`). 33037b752e3dSPatrick Sanan 330460f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCFIELDSPLIT`, `PCFieldSplitGetDetectSaddlePoint()`, `PCFieldSplitSetType()`, `PCFieldSplitSetSchurPre()`, `PC_FIELDSPLIT_SCHUR_PRE_SELF` 33057b752e3dSPatrick Sanan @*/ 3306d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc, PetscBool flg) 3307d71ae5a4SJacob Faibussowitsch { 33087b752e3dSPatrick Sanan PC_FieldSplit *jac = (PC_FieldSplit *)pc->data; 33097b752e3dSPatrick Sanan 33107b752e3dSPatrick Sanan PetscFunctionBegin; 33117b752e3dSPatrick Sanan jac->detect = flg; 33127b752e3dSPatrick Sanan if (jac->detect) { 33139566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR)); 33149566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_SELF, NULL)); 33157b752e3dSPatrick Sanan } 33163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33177b752e3dSPatrick Sanan } 33187b752e3dSPatrick Sanan 33190971522cSBarry Smith /*MC 3320a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 33210b4b7b1cSBarry Smith collections of variables (that may overlap) called fields or splits. Each field often represents a different continuum variable 33220b4b7b1cSBarry Smith represented on a grid, such as velocity, pressure, or temperature. 33230b4b7b1cSBarry Smith In the literature these are sometimes called block preconditioners; but should not be confused with `PCBJACOBI`. 33240b4b7b1cSBarry Smith See [the users manual section on "Solving Block Matrices"](sec_block_matrices) for more details. 33250971522cSBarry Smith 332679416396SBarry Smith Options Database Keys: 3327de37d9f1SPatrick Sanan + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the `%d`'th split 332881540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 3329de37d9f1SPatrick Sanan been supplied explicitly by `-pc_fieldsplit_%d_fields` 333081540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 333180670ca5SBarry Smith when the matrix is not of `MatType` `MATNEST` 3332a51937d4SCarola Kruse . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur,gkb> - type of relaxation or factorization splitting 3333d59693daSPierre Jolivet . -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is `a11`; see `PCFieldSplitSetSchurPre()` 33341d27aa22SBarry Smith . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - set factorization type when using `-pc_fieldsplit_type schur`; 33351d27aa22SBarry Smith see `PCFieldSplitSetSchurFactType()` 333673ff1848SBarry Smith . -pc_fieldsplit_dm_splits <true,false> (default is true) - Whether to use `DMCreateFieldDecomposition()` for splits 3337fb6809a2SPatrick Sanan - -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero diagonal and uses Schur complement with no preconditioner as the solver 333879416396SBarry Smith 3339de37d9f1SPatrick Sanan Options prefixes for inner solvers when using the Schur complement preconditioner are `-fieldsplit_0_` and `-fieldsplit_1_` . 3340de37d9f1SPatrick Sanan The options prefix for the inner solver when using the Golub-Kahan biadiagonalization preconditioner is `-fieldsplit_0_` 334160f59c3bSBarry Smith For all other solvers they are `-fieldsplit_%d_` for the `%d`'th field; use `-fieldsplit_` for all fields. 334260f59c3bSBarry Smith 334322399129SNuno Nobre To set options on the solvers for all blocks, prepend `-fieldsplit_` to all the `PC` 334422399129SNuno Nobre options database keys. For example, `-fieldsplit_pc_type ilu` `-fieldsplit_pc_factor_levels 1`. 334560f59c3bSBarry Smith 334660f59c3bSBarry Smith To set the options on the solvers separate for each block call `PCFieldSplitGetSubKSP()` 334760f59c3bSBarry Smith and set the options directly on the resulting `KSP` object 334860f59c3bSBarry Smith 334960f59c3bSBarry Smith Level: intermediate 33505d4c12cdSJungho Lee 3351c8a0d604SMatthew G Knepley Notes: 335280670ca5SBarry Smith Use `PCFieldSplitSetFields()` to set splits defined by "strided" entries or with a `MATNEST` and `PCFieldSplitSetIS()` 3353f1580f4eSBarry Smith to define a split by an arbitrary collection of entries. 3354d32f9abdSBarry Smith 335573ff1848SBarry Smith If no splits are set, the default is used. If a `DM` is associated with the `PC` and it supports 335680670ca5SBarry Smith `DMCreateFieldDecomposition()`, then that is used for the default. Otherwise if the matrix is not `MATNEST`, the splits are defined by entries strided by bs, 3357de37d9f1SPatrick Sanan beginning at 0 then 1, etc to bs-1. The block size can be set with `PCFieldSplitSetBlockSize()`, 3358d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 3359de37d9f1SPatrick Sanan to `KSPSetOperators()`/`PCSetOperators()`. 3360d32f9abdSBarry Smith 3361de37d9f1SPatrick Sanan For the Schur complement preconditioner if 3362de37d9f1SPatrick Sanan ```{math} 3363de37d9f1SPatrick Sanan J = \left[\begin{array}{cc} A_{00} & A_{01} \\ A_{10} & A_{11} \end{array}\right] 3364de37d9f1SPatrick Sanan ``` 3365e69d4d44SBarry Smith 3366de37d9f1SPatrick Sanan the preconditioner using `full` factorization is logically 3367de37d9f1SPatrick Sanan ```{math} 3368f820a28cSNuno Nobre \left[\begin{array}{cc} I & -\text{ksp}(A_{00}) A_{01} \\ 0 & I \end{array}\right] \left[\begin{array}{cc} \text{ksp}(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] 3369de37d9f1SPatrick Sanan ``` 3370cee94454SNuno Nobre where the action of $\text{ksp}(A_{00})$ is applied using the `KSP` solver with prefix `-fieldsplit_0_`. $S$ is the Schur complement 3371de37d9f1SPatrick Sanan ```{math} 3372223e5b4fSPatrick Sanan S = A_{11} - A_{10} \text{ksp}(A_{00}) A_{01} 3373de37d9f1SPatrick Sanan ``` 33747cbeddf0SNuno Nobre 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` 33757cbeddf0SNuno Nobre was given in providing the SECOND split or 1 if not given). Accordingly, if using `PCFieldSplitGetSubKSP()`, the array of sub-`KSP` contexts will hold two `KSP`s: at its 33767cbeddf0SNuno Nobre 0th index, the `KSP` associated with `-fieldsplit_0_`, and at its 1st index, the `KSP` corresponding to `-fieldsplit_1_`. 3377d82b6cdfSNuno Nobre By default, $A_{11}$ is used to construct a preconditioner for $S$, use `PCFieldSplitSetSchurPre()` for all the possible ways to construct the preconditioner for $S$. 3378de37d9f1SPatrick Sanan 3379de37d9f1SPatrick Sanan The factorization type is set using `-pc_fieldsplit_schur_fact_type <diag, lower, upper, full>`. `full` is shown above, 3380de37d9f1SPatrick Sanan `diag` gives 3381de37d9f1SPatrick Sanan ```{math} 3382f820a28cSNuno Nobre \left[\begin{array}{cc} \text{ksp}(A_{00}) & 0 \\ 0 & -\text{ksp}(S) \end{array}\right] 3383de37d9f1SPatrick Sanan ``` 3384de37d9f1SPatrick 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 3385de37d9f1SPatrick Sanan can be turned off with `PCFieldSplitSetSchurScale()` or by command line `-pc_fieldsplit_schur_scale 1.0`. The `lower` factorization is the inverse of 3386de37d9f1SPatrick Sanan ```{math} 3387de37d9f1SPatrick Sanan \left[\begin{array}{cc} A_{00} & 0 \\ A_{10} & S \end{array}\right] 3388de37d9f1SPatrick Sanan ``` 3389cee94454SNuno Nobre where the inverses of $A_{00}$ and $S$ are applied using `KSP`s. The upper factorization is the inverse of 3390de37d9f1SPatrick Sanan ```{math} 3391de37d9f1SPatrick Sanan \left[\begin{array}{cc} A_{00} & A_{01} \\ 0 & S \end{array}\right] 3392de37d9f1SPatrick Sanan ``` 3393de37d9f1SPatrick Sanan where again the inverses of $A_{00}$ and $S$ are applied using `KSP`s. 3394de37d9f1SPatrick Sanan 3395de37d9f1SPatrick Sanan If only one set of indices (one `IS`) is provided with `PCFieldSplitSetIS()` then the complement of that `IS` 339680670ca5SBarry Smith is used automatically for a second submatrix. 3397edf189efSBarry Smith 3398f1580f4eSBarry Smith The fieldsplit preconditioner cannot currently be used with the `MATBAIJ` or `MATSBAIJ` data formats if the blocksize is larger than 1. 339980670ca5SBarry Smith Generally it should be used with the `MATAIJ` or `MATNEST` `MatType` 3400ff218e97SBarry Smith 340180670ca5SBarry Smith The forms of these preconditioners are closely related, if not identical, to forms derived as "Distributive Iterations", see, 34021d569b8fSBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling {cite}`wesseling2009`. 340380670ca5SBarry Smith One can also use `PCFIELDSPLIT` inside a smoother resulting in "Distributive Smoothers". 34040716a85fSBarry Smith 3405de37d9f1SPatrick Sanan See "A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations" {cite}`elman2008tcp`. 3406a6a584a2SBarry Smith 3407de37d9f1SPatrick Sanan The Constrained Pressure Preconditioner (CPR) can be implemented using `PCCOMPOSITE` with `PCGALERKIN`. CPR first solves an $R A P$ subsystem, updates the 3408de37d9f1SPatrick Sanan residual on all variables (`PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE)`), and then applies a simple ILU like preconditioner on all the variables. 3409a51937d4SCarola Kruse 3410de37d9f1SPatrick Sanan The generalized Golub-Kahan bidiagonalization preconditioner (GKB) can be applied to symmetric $2 \times 2$ block matrices of the shape 3411de37d9f1SPatrick Sanan ```{math} 3412de37d9f1SPatrick Sanan \left[\begin{array}{cc} A_{00} & A_{01} \\ A_{01}' & 0 \end{array}\right] 3413de37d9f1SPatrick Sanan ``` 34141d569b8fSBarry Smith 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}'$. 3415de37d9f1SPatrick 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_`. 3416a51937d4SCarola Kruse 34170b4b7b1cSBarry Smith Some `PCFIELDSPLIT` variants are called physics-based preconditioners, since the preconditioner takes into account the underlying physics of the 34180b4b7b1cSBarry Smith problem. But this nomenclature is not well-defined. 34190b4b7b1cSBarry Smith 342060f59c3bSBarry Smith Developer Note: 342160f59c3bSBarry Smith The Schur complement functionality of `PCFIELDSPLIT` should likely be factored into its own `PC` thus simplifying the implementation of the preconditioners and their 342260f59c3bSBarry Smith user API. 3423f1580f4eSBarry Smith 342460f59c3bSBarry Smith .seealso: [](sec_block_matrices), `PC`, `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCLSC`, 3425db781477SPatrick Sanan `PCFieldSplitGetSubKSP()`, `PCFieldSplitSchurGetSubKSP()`, `PCFieldSplitSetFields()`, 3426db781477SPatrick Sanan `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`, `PCFieldSplitSetSchurFactType()`, 3427db781477SPatrick Sanan `MatSchurComplementSetAinvType()`, `PCFieldSplitSetSchurScale()`, `PCFieldSplitSetDetectSaddlePoint()` 34280971522cSBarry Smith M*/ 34290971522cSBarry Smith 3430d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc) 3431d71ae5a4SJacob Faibussowitsch { 34320971522cSBarry Smith PC_FieldSplit *jac; 34330971522cSBarry Smith 34340971522cSBarry Smith PetscFunctionBegin; 34354dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 34362fa5cd67SKarl Rupp 34370971522cSBarry Smith jac->bs = -1; 34383e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 3439e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 3440c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 3441c096484dSStefano Zampini jac->schurscale = -1.0; 3442fbe7908bSJed Brown jac->dm_splits = PETSC_TRUE; 3443a51937d4SCarola Kruse jac->gkbtol = 1e-5; 3444a51937d4SCarola Kruse jac->gkbdelay = 5; 3445a51937d4SCarola Kruse jac->gkbnu = 1; 3446a51937d4SCarola Kruse jac->gkbmaxit = 100; 344751f519a2SBarry Smith 34480971522cSBarry Smith pc->data = (void *)jac; 34490971522cSBarry Smith 34500971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 3451574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 34520971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 34530971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 34540a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 34550971522cSBarry Smith 34569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSchurGetSubKSP_C", PCFieldSplitSchurGetSubKSP_FieldSplit)); 34579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetFields_C", PCFieldSplitSetFields_FieldSplit)); 34589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetIS_C", PCFieldSplitSetIS_FieldSplit)); 34599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetType_C", PCFieldSplitSetType_FieldSplit)); 34609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitSetBlockSize_C", PCFieldSplitSetBlockSize_FieldSplit)); 34619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFieldSplitRestrictIS_C", PCFieldSplitRestrictIS_FieldSplit)); 34629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_FieldSplit)); 34637ff38633SStefano Zampini 34647ff38633SStefano Zampini /* Initialize function pointers */ 34657ff38633SStefano Zampini PetscCall(PCFieldSplitSetType(pc, jac->type)); 34663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34670971522cSBarry Smith } 3468