1daa6d063SJed Brown 2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 3daa6d063SJed Brown 4daa6d063SJed Brown typedef struct { 5ace3abfcSBarry Smith PetscBool allocated; 6ace3abfcSBarry Smith PetscBool scalediag; 7daa6d063SJed Brown KSP kspL; 8daa6d063SJed Brown Vec scale; 9daa6d063SJed Brown Vec x0,y0,x1; 107e8cb189SBarry Smith Mat L; /* keep a copy to reuse when obtained with L = A10*A01 */ 11daa6d063SJed Brown } PC_LSC; 12daa6d063SJed Brown 13daa6d063SJed Brown static PetscErrorCode PCLSCAllocate_Private(PC pc) 14daa6d063SJed Brown { 15daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 16daa6d063SJed Brown Mat A; 17daa6d063SJed Brown PetscErrorCode ierr; 18daa6d063SJed Brown 19daa6d063SJed Brown PetscFunctionBegin; 20daa6d063SJed Brown if (lsc->allocated) PetscFunctionReturn(0); 21ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&lsc->kspL);CHKERRQ(ierr); 22422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(lsc->kspL,pc->erroriffailure);CHKERRQ(ierr); 237e8cb189SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);CHKERRQ(ierr); 24daa6d063SJed Brown ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr); 25daa6d063SJed Brown ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr); 26daa6d063SJed Brown ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr); 27bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 282a7a6963SBarry Smith ierr = MatCreateVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr); 292a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&lsc->x1,NULL);CHKERRQ(ierr); 30daa6d063SJed Brown if (lsc->scalediag) { 31daa6d063SJed Brown ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr); 32daa6d063SJed Brown } 33daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 34daa6d063SJed Brown PetscFunctionReturn(0); 35daa6d063SJed Brown } 36daa6d063SJed Brown 37daa6d063SJed Brown static PetscErrorCode PCSetUp_LSC(PC pc) 38daa6d063SJed Brown { 39daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 407e8cb189SBarry Smith Mat L,Lp,B,C; 41daa6d063SJed Brown PetscErrorCode ierr; 42daa6d063SJed Brown 43daa6d063SJed Brown PetscFunctionBegin; 44daa6d063SJed Brown ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr); 45093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr); 46093c86ffSJed Brown if (!L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr);} 47daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr); 48093c86ffSJed Brown if (!Lp) {ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr);} 497e8cb189SBarry Smith if (!L) { 50bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr); 517e8cb189SBarry Smith if (!lsc->L) { 527e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 537e8cb189SBarry Smith } else { 547e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 557e8cb189SBarry Smith } 567e8cb189SBarry Smith Lp = L = lsc->L; 577e8cb189SBarry Smith } 58daa6d063SJed Brown if (lsc->scale) { 59daa6d063SJed Brown Mat Ap; 60bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,&Ap,NULL,NULL,NULL);CHKERRQ(ierr); 61daa6d063SJed Brown ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 62daa6d063SJed Brown ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 63daa6d063SJed Brown } 6423ee1639SBarry Smith ierr = KSPSetOperators(lsc->kspL,L,Lp);CHKERRQ(ierr); 65c1146ac7SBarry Smith ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr); 66daa6d063SJed Brown PetscFunctionReturn(0); 67daa6d063SJed Brown } 68daa6d063SJed Brown 69daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 70daa6d063SJed Brown { 71daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 72daa6d063SJed Brown Mat A,B,C; 73daa6d063SJed Brown PetscErrorCode ierr; 74daa6d063SJed Brown 75daa6d063SJed Brown PetscFunctionBegin; 76bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,&B,&C,NULL);CHKERRQ(ierr); 77daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 78*c0decd05SBarry Smith ierr = KSPCheckSolve(lsc->kspL,pc,lsc->x1);CHKERRQ(ierr); 79daa6d063SJed Brown ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 80daa6d063SJed Brown if (lsc->scale) { 81daa6d063SJed Brown ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 82daa6d063SJed Brown } 83daa6d063SJed Brown ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 84daa6d063SJed Brown if (lsc->scale) { 85daa6d063SJed Brown ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 86daa6d063SJed Brown } 87daa6d063SJed Brown ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 88daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 89*c0decd05SBarry Smith ierr = KSPCheckSolve(lsc->kspL,pc,y);CHKERRQ(ierr); 90daa6d063SJed Brown PetscFunctionReturn(0); 91daa6d063SJed Brown } 92daa6d063SJed Brown 93a06653b4SBarry Smith static PetscErrorCode PCReset_LSC(PC pc) 94daa6d063SJed Brown { 95daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 96daa6d063SJed Brown PetscErrorCode ierr; 97daa6d063SJed Brown 98daa6d063SJed Brown PetscFunctionBegin; 996bf464f9SBarry Smith ierr = VecDestroy(&lsc->x0);CHKERRQ(ierr); 1006bf464f9SBarry Smith ierr = VecDestroy(&lsc->y0);CHKERRQ(ierr); 1016bf464f9SBarry Smith ierr = VecDestroy(&lsc->x1);CHKERRQ(ierr); 1026bf464f9SBarry Smith ierr = VecDestroy(&lsc->scale);CHKERRQ(ierr); 1036bf464f9SBarry Smith ierr = KSPDestroy(&lsc->kspL);CHKERRQ(ierr); 1046bf464f9SBarry Smith ierr = MatDestroy(&lsc->L);CHKERRQ(ierr); 105a06653b4SBarry Smith PetscFunctionReturn(0); 106a06653b4SBarry Smith } 107a06653b4SBarry Smith 108a06653b4SBarry Smith static PetscErrorCode PCDestroy_LSC(PC pc) 109a06653b4SBarry Smith { 110a06653b4SBarry Smith PetscErrorCode ierr; 111a06653b4SBarry Smith 112a06653b4SBarry Smith PetscFunctionBegin; 113a06653b4SBarry Smith ierr = PCReset_LSC(pc);CHKERRQ(ierr); 114c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 115daa6d063SJed Brown PetscFunctionReturn(0); 116daa6d063SJed Brown } 117daa6d063SJed Brown 1184416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LSC(PetscOptionItems *PetscOptionsObject,PC pc) 119daa6d063SJed Brown { 120daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 121daa6d063SJed Brown PetscErrorCode ierr; 122daa6d063SJed Brown 123daa6d063SJed Brown PetscFunctionBegin; 124e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"LSC options");CHKERRQ(ierr); 125e55864a3SBarry Smith { 1260298fd71SBarry Smith ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,NULL);CHKERRQ(ierr); 1271a1499c8SBarry Smith } 128daa6d063SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 129daa6d063SJed Brown PetscFunctionReturn(0); 130daa6d063SJed Brown } 131daa6d063SJed Brown 1327e8cb189SBarry Smith static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 1337e8cb189SBarry Smith { 1347e8cb189SBarry Smith PC_LSC *jac = (PC_LSC*)pc->data; 1357e8cb189SBarry Smith PetscErrorCode ierr; 1367e8cb189SBarry Smith PetscBool iascii; 1377e8cb189SBarry Smith 1387e8cb189SBarry Smith PetscFunctionBegin; 139251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1407e8cb189SBarry Smith if (iascii) { 1417e8cb189SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1427e8cb189SBarry Smith ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 1437e8cb189SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 14411aeaf0aSBarry Smith } 1457e8cb189SBarry Smith PetscFunctionReturn(0); 1467e8cb189SBarry Smith } 1477e8cb189SBarry Smith 148daa6d063SJed Brown /*MC 149daa6d063SJed Brown PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 150daa6d063SJed Brown 151daa6d063SJed Brown Options Database Key: 152daa6d063SJed Brown . -pc_lsc_scale_diag - Use the diagonal of A for scaling 153daa6d063SJed Brown 154daa6d063SJed Brown Level: intermediate 155daa6d063SJed Brown 156daa6d063SJed Brown Notes: 157daa6d063SJed Brown This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 158daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 159daa6d063SJed Brown 160daa6d063SJed Brown .vb 1617e8cb189SBarry Smith S = A11 - A10 inv(A00) A01 162daa6d063SJed Brown .ve 163daa6d063SJed Brown 1647e8cb189SBarry Smith PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 165daa6d063SJed Brown inv(S) is given by 166daa6d063SJed Brown 167daa6d063SJed Brown .vb 1687e8cb189SBarry Smith inv(A10 A01) A10 A00 A01 inv(A10 A01) 169daa6d063SJed Brown .ve 170daa6d063SJed Brown 1717e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 1727e8cb189SBarry Smith usually more efficient anyway). In the case of incompressible flow, A10 A10 is a Laplacian, call it L. The current 173daa6d063SJed Brown interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 174daa6d063SJed Brown 17562279869SBarry Smith If you had called KSPSetOperators(ksp,S,Sp), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 176daa6d063SJed Brown like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 177daa6d063SJed Brown For example, you might have setup code like this 178daa6d063SJed Brown 179daa6d063SJed Brown .vb 180daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 181daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 182daa6d063SJed Brown .ve 183daa6d063SJed Brown 184daa6d063SJed Brown And then your Jacobian assembly would look like 185daa6d063SJed Brown 186daa6d063SJed Brown .vb 187daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 188daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 189daa6d063SJed Brown if (L) { assembly L } 190daa6d063SJed Brown if (Lp) { assemble Lp } 191daa6d063SJed Brown .ve 192daa6d063SJed Brown 193daa6d063SJed Brown With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 194daa6d063SJed Brown 195daa6d063SJed Brown .vb 196daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 197daa6d063SJed Brown .ve 198daa6d063SJed Brown 199daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 200daa6d063SJed Brown 201a674a499SJed Brown References: 20296a0c994SBarry Smith + 1. - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006. 20396a0c994SBarry Smith - 2. - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001. 204a674a499SJed Brown 205daa6d063SJed Brown Concepts: physics based preconditioners, block preconditioners 206daa6d063SJed Brown 207daa6d063SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 20829f8a81cSJed Brown PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(), 209daa6d063SJed Brown MatCreateSchurComplement() 210daa6d063SJed Brown M*/ 211daa6d063SJed Brown 2128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc) 213daa6d063SJed Brown { 214daa6d063SJed Brown PC_LSC *lsc; 215daa6d063SJed Brown PetscErrorCode ierr; 216daa6d063SJed Brown 217daa6d063SJed Brown PetscFunctionBegin; 218b00a9115SJed Brown ierr = PetscNewLog(pc,&lsc);CHKERRQ(ierr); 219daa6d063SJed Brown pc->data = (void*)lsc; 220daa6d063SJed Brown 221daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 222daa6d063SJed Brown pc->ops->applytranspose = 0; 223daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 224a06653b4SBarry Smith pc->ops->reset = PCReset_LSC; 225daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 226daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 2277e8cb189SBarry Smith pc->ops->view = PCView_LSC; 228daa6d063SJed Brown pc->ops->applyrichardson = 0; 229daa6d063SJed Brown PetscFunctionReturn(0); 230daa6d063SJed Brown } 231