1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 2daa6d063SJed Brown 3daa6d063SJed Brown typedef struct { 49fd865deSAlexander PetscBool allocated, commute, scalediag; 59fd865deSAlexander KSP kspL, kspMass; 69fd865deSAlexander Vec Avec0, Avec1, Svec0, scale; 79fd865deSAlexander Mat L; 8daa6d063SJed Brown } PC_LSC; 9daa6d063SJed Brown 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCLSCAllocate_Private(PC pc) 11d71ae5a4SJacob Faibussowitsch { 12daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 13daa6d063SJed Brown Mat A; 14daa6d063SJed Brown 15daa6d063SJed Brown PetscFunctionBegin; 163ba16761SJacob Faibussowitsch if (lsc->allocated) PetscFunctionReturn(PETSC_SUCCESS); 179566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspL)); 183821be0aSBarry Smith PetscCall(KSPSetNestLevel(lsc->kspL, pc->kspnestlevel)); 199566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(lsc->kspL, pc->erroriffailure)); 209566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)lsc->kspL, (PetscObject)pc, 1)); 219566063dSJacob Faibussowitsch PetscCall(KSPSetType(lsc->kspL, KSPPREONLY)); 229566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(lsc->kspL, ((PetscObject)pc)->prefix)); 239566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(lsc->kspL, "lsc_")); 249566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, NULL, NULL, NULL)); 259fd865deSAlexander PetscCall(MatCreateVecs(A, &lsc->Avec0, &lsc->Avec1)); 269fd865deSAlexander PetscCall(MatCreateVecs(pc->pmat, &lsc->Svec0, NULL)); 279fd865deSAlexander if (lsc->scalediag) PetscCall(VecDuplicate(lsc->Avec0, &lsc->scale)); 289fd865deSAlexander 299fd865deSAlexander if (lsc->commute) { 309fd865deSAlexander PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspMass)); 319fd865deSAlexander PetscCall(KSPSetErrorIfNotConverged(lsc->kspMass, pc->erroriffailure)); 329fd865deSAlexander PetscCall(PetscObjectIncrementTabLevel((PetscObject)lsc->kspMass, (PetscObject)pc, 1)); 339fd865deSAlexander PetscCall(KSPSetType(lsc->kspMass, KSPPREONLY)); 349fd865deSAlexander PetscCall(KSPSetOptionsPrefix(lsc->kspMass, ((PetscObject)pc)->prefix)); 359fd865deSAlexander PetscCall(KSPAppendOptionsPrefix(lsc->kspMass, "lsc_mass_")); 369fd865deSAlexander } else lsc->kspMass = NULL; 379fd865deSAlexander 38daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40daa6d063SJed Brown } 41daa6d063SJed Brown 42d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LSC(PC pc) 43d71ae5a4SJacob Faibussowitsch { 44daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 459fd865deSAlexander Mat L, Lp, Qscale; 46daa6d063SJed Brown 47daa6d063SJed Brown PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(PCLSCAllocate_Private(pc)); 499fd865deSAlexander 509fd865deSAlexander /* Query for L operator */ 519566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_L", (PetscObject *)&L)); 529566063dSJacob Faibussowitsch if (!L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_L", (PetscObject *)&L)); 539566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Lp", (PetscObject *)&Lp)); 549566063dSJacob Faibussowitsch if (!Lp) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Lp", (PetscObject *)&Lp)); 559fd865deSAlexander 569fd865deSAlexander /* Query for mass operator */ 579fd865deSAlexander PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Qscale", (PetscObject *)&Qscale)); 589fd865deSAlexander if (!Qscale) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Qscale", (PetscObject *)&Qscale)); 599fd865deSAlexander 609fd865deSAlexander if (lsc->commute) { 619fd865deSAlexander PetscCheck(L || Lp, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "The user must provide an L operator for LSC preconditioning when commuting"); 629fd865deSAlexander if (!L && Lp) L = Lp; 639fd865deSAlexander else if (L && !Lp) Lp = L; 649fd865deSAlexander 659fd865deSAlexander PetscCheck(Qscale, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "The user must provide a Qscale matrix for LSC preconditioning when commuting"); 667e8cb189SBarry Smith } else { 679fd865deSAlexander if (lsc->scale) { 689fd865deSAlexander if (!Qscale) PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, &Qscale, NULL, NULL, NULL)); 699fd865deSAlexander PetscCall(MatGetDiagonal(Qscale, lsc->scale)); 709fd865deSAlexander PetscCall(VecReciprocal(lsc->scale)); 719fd865deSAlexander } 729fd865deSAlexander if (!L) { 739fd865deSAlexander Mat B, C; 749fd865deSAlexander PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, NULL, &B, &C, NULL)); 759fd865deSAlexander if (lsc->scale) { 769fd865deSAlexander Mat CAdiaginv; 779fd865deSAlexander PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &CAdiaginv)); 789fd865deSAlexander PetscCall(MatDiagonalScale(CAdiaginv, NULL, lsc->scale)); 799fd865deSAlexander if (!lsc->L) PetscCall(MatMatMult(CAdiaginv, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &lsc->L)); 809fd865deSAlexander else PetscCall(MatMatMult(CAdiaginv, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &lsc->L)); 819fd865deSAlexander PetscCall(MatDestroy(&CAdiaginv)); 829fd865deSAlexander } else { 839fd865deSAlexander if (!lsc->L) { 849fd865deSAlexander PetscCall(MatProductCreate(C, B, NULL, &lsc->L)); 859fd865deSAlexander PetscCall(MatProductSetType(lsc->L, MATPRODUCT_AB)); 869fd865deSAlexander PetscCall(MatProductSetFromOptions(lsc->L)); 879fd865deSAlexander PetscCall(MatProductSymbolic(lsc->L)); 889fd865deSAlexander } 899fd865deSAlexander PetscCall(MatProductNumeric(lsc->L)); 907e8cb189SBarry Smith } 917e8cb189SBarry Smith Lp = L = lsc->L; 927e8cb189SBarry Smith } 93daa6d063SJed Brown } 949fd865deSAlexander 959566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(lsc->kspL, L, Lp)); 969566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(lsc->kspL)); 979fd865deSAlexander if (lsc->commute) { 989fd865deSAlexander PetscCall(KSPSetOperators(lsc->kspMass, Qscale, Qscale)); 999fd865deSAlexander PetscCall(KSPSetFromOptions(lsc->kspMass)); 1009fd865deSAlexander } 1013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102daa6d063SJed Brown } 103daa6d063SJed Brown 104d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LSC(PC pc, Vec x, Vec y) 105d71ae5a4SJacob Faibussowitsch { 106daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 107daa6d063SJed Brown Mat A, B, C; 108daa6d063SJed Brown 109daa6d063SJed Brown PetscFunctionBegin; 1109566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, &B, &C, NULL)); 1119fd865deSAlexander if (lsc->commute) { 1129fd865deSAlexander PetscCall(KSPSolve(lsc->kspMass, x, lsc->Svec0)); 1139fd865deSAlexander PetscCall(KSPCheckSolve(lsc->kspMass, pc, lsc->Svec0)); 1149fd865deSAlexander PetscCall(MatMult(B, lsc->Svec0, lsc->Avec0)); 1159fd865deSAlexander PetscCall(KSPSolve(lsc->kspL, lsc->Avec0, lsc->Avec1)); 1169fd865deSAlexander PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->Avec1)); 1179fd865deSAlexander PetscCall(MatMult(A, lsc->Avec1, lsc->Avec0)); 1189fd865deSAlexander PetscCall(KSPSolve(lsc->kspL, lsc->Avec0, lsc->Avec1)); 1199fd865deSAlexander PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->Avec1)); 1209fd865deSAlexander PetscCall(MatMult(C, lsc->Avec1, lsc->Svec0)); 1219fd865deSAlexander PetscCall(KSPSolve(lsc->kspMass, lsc->Svec0, y)); 1229fd865deSAlexander PetscCall(KSPCheckSolve(lsc->kspMass, pc, y)); 1239fd865deSAlexander } else { 1249fd865deSAlexander PetscCall(KSPSolve(lsc->kspL, x, lsc->Svec0)); 1259fd865deSAlexander PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->Svec0)); 1269fd865deSAlexander PetscCall(MatMult(B, lsc->Svec0, lsc->Avec0)); 1279fd865deSAlexander if (lsc->scale) PetscCall(VecPointwiseMult(lsc->Avec0, lsc->Avec0, lsc->scale)); 1289fd865deSAlexander PetscCall(MatMult(A, lsc->Avec0, lsc->Avec1)); 1299fd865deSAlexander if (lsc->scale) PetscCall(VecPointwiseMult(lsc->Avec1, lsc->Avec1, lsc->scale)); 1309fd865deSAlexander PetscCall(MatMult(C, lsc->Avec1, lsc->Svec0)); 1319fd865deSAlexander PetscCall(KSPSolve(lsc->kspL, lsc->Svec0, y)); 1329566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(lsc->kspL, pc, y)); 1339fd865deSAlexander } 1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135daa6d063SJed Brown } 136daa6d063SJed Brown 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LSC(PC pc) 138d71ae5a4SJacob Faibussowitsch { 139daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 140daa6d063SJed Brown 141daa6d063SJed Brown PetscFunctionBegin; 1429fd865deSAlexander PetscCall(VecDestroy(&lsc->Avec0)); 1439fd865deSAlexander PetscCall(VecDestroy(&lsc->Avec1)); 1449fd865deSAlexander PetscCall(VecDestroy(&lsc->Svec0)); 1459566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&lsc->kspL)); 1469fd865deSAlexander if (lsc->commute) PetscCall(KSPDestroy(&lsc->kspMass)); 1479fd865deSAlexander if (lsc->L) PetscCall(MatDestroy(&lsc->L)); 1489fd865deSAlexander if (lsc->scale) PetscCall(VecDestroy(&lsc->scale)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150a06653b4SBarry Smith } 151a06653b4SBarry Smith 152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LSC(PC pc) 153d71ae5a4SJacob Faibussowitsch { 154a06653b4SBarry Smith PetscFunctionBegin; 1559566063dSJacob Faibussowitsch PetscCall(PCReset_LSC(pc)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158daa6d063SJed Brown } 159daa6d063SJed Brown 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_LSC(PC pc, PetscOptionItems *PetscOptionsObject) 161d71ae5a4SJacob Faibussowitsch { 162daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 163daa6d063SJed Brown 164daa6d063SJed Brown PetscFunctionBegin; 165d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "LSC options"); 166d71ae5a4SJacob Faibussowitsch { 1679fd865deSAlexander PetscCall(PetscOptionsBool("-pc_lsc_commute", "Whether to commute the LSC preconditioner in the style of Olshanskii", "None", lsc->commute, &lsc->commute, NULL)); 1689fd865deSAlexander PetscCall(PetscOptionsBool("-pc_lsc_scale_diag", "Whether to scale BBt products. Will use the inverse of the diagonal of Qscale or A if the former is not provided.", "None", lsc->scalediag, &lsc->scalediag, NULL)); 1699fd865deSAlexander PetscCheck(!lsc->scalediag || !lsc->commute, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Diagonal-based scaling is not used when doing a commuted LSC. Either do not ask for diagonal-based scaling or use non-commuted LSC in the original style of Elman"); 170d71ae5a4SJacob Faibussowitsch } 171d0609cedSBarry Smith PetscOptionsHeadEnd(); 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173daa6d063SJed Brown } 174daa6d063SJed Brown 175d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_LSC(PC pc, PetscViewer viewer) 176d71ae5a4SJacob Faibussowitsch { 1777e8cb189SBarry Smith PC_LSC *jac = (PC_LSC *)pc->data; 1787e8cb189SBarry Smith PetscBool iascii; 1797e8cb189SBarry Smith 1807e8cb189SBarry Smith PetscFunctionBegin; 1819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1827e8cb189SBarry Smith if (iascii) { 1839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 184ca0c957dSBarry Smith if (jac->kspL) { 1859566063dSJacob Faibussowitsch PetscCall(KSPView(jac->kspL, viewer)); 186ca0c957dSBarry Smith } else { 1879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC KSP object not yet created, hence cannot display")); 188ca0c957dSBarry Smith } 1899fd865deSAlexander if (jac->commute) { 1909fd865deSAlexander if (jac->kspMass) { 1919fd865deSAlexander PetscCall(KSPView(jac->kspMass, viewer)); 1929fd865deSAlexander } else { 1939fd865deSAlexander PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC Mass KSP object not yet created, hence cannot display")); 1949fd865deSAlexander } 1959fd865deSAlexander } 1969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 19711aeaf0aSBarry Smith } 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1997e8cb189SBarry Smith } 2007e8cb189SBarry Smith 201daa6d063SJed Brown /*MC 202*1d27aa22SBarry Smith PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators {cite}`elmanhowleshadidshuttleworthtuminaro2006` {cite}`silvester2001efficient` 203daa6d063SJed Brown 204daa6d063SJed Brown Options Database Key: 205c3ff4b9fSAlex Lindsay + -pc_lsc_commute - Whether to commute the LSC preconditioner in the style of Olshanskii 206*1d27aa22SBarry Smith - -pc_lsc_scale_diag - Whether to scale $BB^T$ products. Will use the inverse of the diagonal of Qscale or A if the former is not provided 207daa6d063SJed Brown 208daa6d063SJed Brown Level: intermediate 209daa6d063SJed Brown 210daa6d063SJed Brown Notes: 211f1580f4eSBarry Smith This preconditioner will normally be used with `PCFIELDSPLIT` to precondition the Schur complement, but 212daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 213daa6d063SJed Brown 214*1d27aa22SBarry Smith $$ 215*1d27aa22SBarry Smith S = A11 - A10 A00^{-1} A01 216*1d27aa22SBarry Smith $$ 217daa6d063SJed Brown 218f1580f4eSBarry Smith `PCLSC` currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 219daa6d063SJed Brown inv(S) is given by 220daa6d063SJed Brown 221*1d27aa22SBarry Smith $$ 222*1d27aa22SBarry Smith (A10 A01)^{-1} A10 A00 A01 (A10 A01)^{-1} 223*1d27aa22SBarry Smith $$ 224daa6d063SJed Brown 2257e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 2261c6f1693SPatrick Sanan usually more efficient anyway). In the case of incompressible flow, A10 A01 is a Laplacian; call it L. The current 227*1d27aa22SBarry Smith interface is to compose L and a preconditioning matrix Lp on the preconditioning matrix. 228daa6d063SJed Brown 229f1580f4eSBarry Smith If you had called `KSPSetOperators`(ksp,S,Sp), S should have type `MATSCHURCOMPLEMENT` and Sp can be any type you 230f1580f4eSBarry Smith like (`PCLSC` doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 231daa6d063SJed Brown For example, you might have setup code like this 232daa6d063SJed Brown 233daa6d063SJed Brown .vb 234daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 235daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 236daa6d063SJed Brown .ve 237daa6d063SJed Brown 238daa6d063SJed Brown And then your Jacobian assembly would look like 239daa6d063SJed Brown 240daa6d063SJed Brown .vb 241daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 242daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 243daa6d063SJed Brown if (L) { assembly L } 244daa6d063SJed Brown if (Lp) { assemble Lp } 245daa6d063SJed Brown .ve 246daa6d063SJed Brown 247*1d27aa22SBarry Smith With this, you should be able to choose LSC preconditioning, using e.g. the `PCML` algebraic multigrid to solve with L 248daa6d063SJed Brown .vb 249daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 250daa6d063SJed Brown .ve 251daa6d063SJed Brown 252daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 253daa6d063SJed Brown 254562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `Block_Preconditioners`, `PCFIELDSPLIT`, 255db781477SPatrick Sanan `PCFieldSplitGetSubKSP()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`, 256f1580f4eSBarry Smith `MatCreateSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetSubMatrices()`, `MatSchurComplementUpdateSubMatrices()`, 257f1580f4eSBarry Smith `MatSchurComplementSetAinvType()`, `MatGetSchurComplement()` 258daa6d063SJed Brown M*/ 259daa6d063SJed Brown 260d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc) 261d71ae5a4SJacob Faibussowitsch { 262daa6d063SJed Brown PC_LSC *lsc; 263daa6d063SJed Brown 264daa6d063SJed Brown PetscFunctionBegin; 2654dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lsc)); 266daa6d063SJed Brown pc->data = (void *)lsc; 267daa6d063SJed Brown 268daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 2690a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 270daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 271a06653b4SBarry Smith pc->ops->reset = PCReset_LSC; 272daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 273daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 2747e8cb189SBarry Smith pc->ops->view = PCView_LSC; 2750a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277daa6d063SJed Brown } 278