1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 2daa6d063SJed Brown 3daa6d063SJed Brown typedef struct { 4ace3abfcSBarry Smith PetscBool allocated; 5ace3abfcSBarry Smith PetscBool scalediag; 6daa6d063SJed Brown KSP kspL; 7daa6d063SJed Brown Vec scale; 8daa6d063SJed Brown Vec x0, y0, x1; 97e8cb189SBarry Smith Mat L; /* keep a copy to reuse when obtained with L = A10*A01 */ 10daa6d063SJed Brown } PC_LSC; 11daa6d063SJed Brown 129371c9d4SSatish Balay static PetscErrorCode PCLSCAllocate_Private(PC pc) { 13daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 14daa6d063SJed Brown Mat A; 15daa6d063SJed Brown 16daa6d063SJed Brown PetscFunctionBegin; 17daa6d063SJed Brown if (lsc->allocated) PetscFunctionReturn(0); 189566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspL)); 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)); 259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lsc->x0, &lsc->y0)); 269566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &lsc->x1, NULL)); 2748a46eb9SPierre Jolivet if (lsc->scalediag) PetscCall(VecDuplicate(lsc->x0, &lsc->scale)); 28daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 29daa6d063SJed Brown PetscFunctionReturn(0); 30daa6d063SJed Brown } 31daa6d063SJed Brown 329371c9d4SSatish Balay static PetscErrorCode PCSetUp_LSC(PC pc) { 33daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 347e8cb189SBarry Smith Mat L, Lp, B, C; 35daa6d063SJed Brown 36daa6d063SJed Brown PetscFunctionBegin; 379566063dSJacob Faibussowitsch PetscCall(PCLSCAllocate_Private(pc)); 389566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_L", (PetscObject *)&L)); 399566063dSJacob Faibussowitsch if (!L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_L", (PetscObject *)&L)); 409566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Lp", (PetscObject *)&Lp)); 419566063dSJacob Faibussowitsch if (!Lp) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Lp", (PetscObject *)&Lp)); 427e8cb189SBarry Smith if (!L) { 439566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, NULL, &B, &C, NULL)); 447e8cb189SBarry Smith if (!lsc->L) { 459566063dSJacob Faibussowitsch PetscCall(MatMatMult(C, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &lsc->L)); 467e8cb189SBarry Smith } else { 479566063dSJacob Faibussowitsch PetscCall(MatMatMult(C, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &lsc->L)); 487e8cb189SBarry Smith } 497e8cb189SBarry Smith Lp = L = lsc->L; 507e8cb189SBarry Smith } 51daa6d063SJed Brown if (lsc->scale) { 52daa6d063SJed Brown Mat Ap; 539566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, &Ap, NULL, NULL, NULL)); 549566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Ap, lsc->scale)); /* Should be the mass matrix, but we don't have plumbing for that yet */ 559566063dSJacob Faibussowitsch PetscCall(VecReciprocal(lsc->scale)); 56daa6d063SJed Brown } 579566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(lsc->kspL, L, Lp)); 589566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(lsc->kspL)); 59daa6d063SJed Brown PetscFunctionReturn(0); 60daa6d063SJed Brown } 61daa6d063SJed Brown 629371c9d4SSatish Balay static PetscErrorCode PCApply_LSC(PC pc, Vec x, Vec y) { 63daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 64daa6d063SJed Brown Mat A, B, C; 65daa6d063SJed Brown 66daa6d063SJed Brown PetscFunctionBegin; 679566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, &B, &C, NULL)); 689566063dSJacob Faibussowitsch PetscCall(KSPSolve(lsc->kspL, x, lsc->x1)); 699566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->x1)); 709566063dSJacob Faibussowitsch PetscCall(MatMult(B, lsc->x1, lsc->x0)); 711baa6e33SBarry Smith if (lsc->scale) PetscCall(VecPointwiseMult(lsc->x0, lsc->x0, lsc->scale)); 729566063dSJacob Faibussowitsch PetscCall(MatMult(A, lsc->x0, lsc->y0)); 731baa6e33SBarry Smith if (lsc->scale) PetscCall(VecPointwiseMult(lsc->y0, lsc->y0, lsc->scale)); 749566063dSJacob Faibussowitsch PetscCall(MatMult(C, lsc->y0, lsc->x1)); 759566063dSJacob Faibussowitsch PetscCall(KSPSolve(lsc->kspL, lsc->x1, y)); 769566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(lsc->kspL, pc, y)); 77daa6d063SJed Brown PetscFunctionReturn(0); 78daa6d063SJed Brown } 79daa6d063SJed Brown 809371c9d4SSatish Balay static PetscErrorCode PCReset_LSC(PC pc) { 81daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 82daa6d063SJed Brown 83daa6d063SJed Brown PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lsc->x0)); 859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lsc->y0)); 869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lsc->x1)); 879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lsc->scale)); 889566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&lsc->kspL)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lsc->L)); 90a06653b4SBarry Smith PetscFunctionReturn(0); 91a06653b4SBarry Smith } 92a06653b4SBarry Smith 939371c9d4SSatish Balay static PetscErrorCode PCDestroy_LSC(PC pc) { 94a06653b4SBarry Smith PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(PCReset_LSC(pc)); 969566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 97daa6d063SJed Brown PetscFunctionReturn(0); 98daa6d063SJed Brown } 99daa6d063SJed Brown 1009371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_LSC(PC pc, PetscOptionItems *PetscOptionsObject) { 101daa6d063SJed Brown PC_LSC *lsc = (PC_LSC *)pc->data; 102daa6d063SJed Brown 103daa6d063SJed Brown PetscFunctionBegin; 104d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "LSC options"); 1059371c9d4SSatish Balay { PetscCall(PetscOptionsBool("-pc_lsc_scale_diag", "Use diagonal of velocity block (A) for scaling", "None", lsc->scalediag, &lsc->scalediag, NULL)); } 106d0609cedSBarry Smith PetscOptionsHeadEnd(); 107daa6d063SJed Brown PetscFunctionReturn(0); 108daa6d063SJed Brown } 109daa6d063SJed Brown 1109371c9d4SSatish Balay static PetscErrorCode PCView_LSC(PC pc, PetscViewer viewer) { 1117e8cb189SBarry Smith PC_LSC *jac = (PC_LSC *)pc->data; 1127e8cb189SBarry Smith PetscBool iascii; 1137e8cb189SBarry Smith 1147e8cb189SBarry Smith PetscFunctionBegin; 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1167e8cb189SBarry Smith if (iascii) { 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 118ca0c957dSBarry Smith if (jac->kspL) { 1199566063dSJacob Faibussowitsch PetscCall(KSPView(jac->kspL, viewer)); 120ca0c957dSBarry Smith } else { 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC KSP object not yet created, hence cannot display")); 122ca0c957dSBarry Smith } 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 12411aeaf0aSBarry Smith } 1257e8cb189SBarry Smith PetscFunctionReturn(0); 1267e8cb189SBarry Smith } 1277e8cb189SBarry Smith 128daa6d063SJed Brown /*MC 129daa6d063SJed Brown PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 130daa6d063SJed Brown 131daa6d063SJed Brown Options Database Key: 132daa6d063SJed Brown . -pc_lsc_scale_diag - Use the diagonal of A for scaling 133daa6d063SJed Brown 134daa6d063SJed Brown Level: intermediate 135daa6d063SJed Brown 136daa6d063SJed Brown Notes: 137*f1580f4eSBarry Smith This preconditioner will normally be used with `PCFIELDSPLIT` to precondition the Schur complement, but 138daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 139daa6d063SJed Brown 140daa6d063SJed Brown .vb 1417e8cb189SBarry Smith S = A11 - A10 inv(A00) A01 142daa6d063SJed Brown .ve 143daa6d063SJed Brown 144*f1580f4eSBarry Smith `PCLSC` currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 145daa6d063SJed Brown inv(S) is given by 146daa6d063SJed Brown 147daa6d063SJed Brown .vb 1487e8cb189SBarry Smith inv(A10 A01) A10 A00 A01 inv(A10 A01) 149daa6d063SJed Brown .ve 150daa6d063SJed Brown 1517e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 1521c6f1693SPatrick Sanan usually more efficient anyway). In the case of incompressible flow, A10 A01 is a Laplacian; call it L. The current 153daa6d063SJed Brown interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 154daa6d063SJed Brown 155*f1580f4eSBarry Smith If you had called `KSPSetOperators`(ksp,S,Sp), S should have type `MATSCHURCOMPLEMENT` and Sp can be any type you 156*f1580f4eSBarry Smith like (`PCLSC` doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 157daa6d063SJed Brown For example, you might have setup code like this 158daa6d063SJed Brown 159daa6d063SJed Brown .vb 160daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 161daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 162daa6d063SJed Brown .ve 163daa6d063SJed Brown 164daa6d063SJed Brown And then your Jacobian assembly would look like 165daa6d063SJed Brown 166daa6d063SJed Brown .vb 167daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 168daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 169daa6d063SJed Brown if (L) { assembly L } 170daa6d063SJed Brown if (Lp) { assemble Lp } 171daa6d063SJed Brown .ve 172daa6d063SJed Brown 173daa6d063SJed Brown With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 174daa6d063SJed Brown 175daa6d063SJed Brown .vb 176daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 177daa6d063SJed Brown .ve 178daa6d063SJed Brown 179daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 180daa6d063SJed Brown 181a674a499SJed Brown References: 182606c0280SSatish Balay + * - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006. 183606c0280SSatish Balay - * - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001. 184a674a499SJed Brown 185db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `Block_Preconditioners`, `PCFIELDSPLIT`, 186db781477SPatrick Sanan `PCFieldSplitGetSubKSP()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`, 187*f1580f4eSBarry Smith `MatCreateSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetSubMatrices()`, `MatSchurComplementUpdateSubMatrices()`, 188*f1580f4eSBarry Smith `MatSchurComplementSetAinvType()`, `MatGetSchurComplement()` 189daa6d063SJed Brown M*/ 190daa6d063SJed Brown 1919371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc) { 192daa6d063SJed Brown PC_LSC *lsc; 193daa6d063SJed Brown 194daa6d063SJed Brown PetscFunctionBegin; 1959566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &lsc)); 196daa6d063SJed Brown pc->data = (void *)lsc; 197daa6d063SJed Brown 198daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 1990a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 200daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 201a06653b4SBarry Smith pc->ops->reset = PCReset_LSC; 202daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 203daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 2047e8cb189SBarry Smith pc->ops->view = PCView_LSC; 2050a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 206daa6d063SJed Brown PetscFunctionReturn(0); 207daa6d063SJed Brown } 208