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 12daa6d063SJed Brown static PetscErrorCode PCLSCAllocate_Private(PC pc) 13daa6d063SJed Brown { 14daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 15daa6d063SJed Brown Mat A; 16daa6d063SJed Brown PetscErrorCode ierr; 17daa6d063SJed Brown 18daa6d063SJed Brown PetscFunctionBegin; 19daa6d063SJed Brown if (lsc->allocated) PetscFunctionReturn(0); 20ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&lsc->kspL);CHKERRQ(ierr); 21422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(lsc->kspL,pc->erroriffailure);CHKERRQ(ierr); 227e8cb189SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);CHKERRQ(ierr); 23daa6d063SJed Brown ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr); 24daa6d063SJed Brown ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr); 25daa6d063SJed Brown ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr); 26bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 272a7a6963SBarry Smith ierr = MatCreateVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr); 282a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&lsc->x1,NULL);CHKERRQ(ierr); 29daa6d063SJed Brown if (lsc->scalediag) { 30daa6d063SJed Brown ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr); 31daa6d063SJed Brown } 32daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 33daa6d063SJed Brown PetscFunctionReturn(0); 34daa6d063SJed Brown } 35daa6d063SJed Brown 36daa6d063SJed Brown static PetscErrorCode PCSetUp_LSC(PC pc) 37daa6d063SJed Brown { 38daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 397e8cb189SBarry Smith Mat L,Lp,B,C; 40daa6d063SJed Brown PetscErrorCode ierr; 41daa6d063SJed Brown 42daa6d063SJed Brown PetscFunctionBegin; 43daa6d063SJed Brown ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr); 44093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr); 45093c86ffSJed Brown if (!L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr);} 46daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr); 47093c86ffSJed Brown if (!Lp) {ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr);} 487e8cb189SBarry Smith if (!L) { 49bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr); 507e8cb189SBarry Smith if (!lsc->L) { 517e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 527e8cb189SBarry Smith } else { 537e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 547e8cb189SBarry Smith } 557e8cb189SBarry Smith Lp = L = lsc->L; 567e8cb189SBarry Smith } 57daa6d063SJed Brown if (lsc->scale) { 58daa6d063SJed Brown Mat Ap; 59bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,&Ap,NULL,NULL,NULL);CHKERRQ(ierr); 60daa6d063SJed Brown ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 61daa6d063SJed Brown ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 62daa6d063SJed Brown } 6323ee1639SBarry Smith ierr = KSPSetOperators(lsc->kspL,L,Lp);CHKERRQ(ierr); 64c1146ac7SBarry Smith ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr); 65daa6d063SJed Brown PetscFunctionReturn(0); 66daa6d063SJed Brown } 67daa6d063SJed Brown 68daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 69daa6d063SJed Brown { 70daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 71daa6d063SJed Brown Mat A,B,C; 72daa6d063SJed Brown PetscErrorCode ierr; 73daa6d063SJed Brown 74daa6d063SJed Brown PetscFunctionBegin; 75bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,&B,&C,NULL);CHKERRQ(ierr); 76daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 77c0decd05SBarry Smith ierr = KSPCheckSolve(lsc->kspL,pc,lsc->x1);CHKERRQ(ierr); 78daa6d063SJed Brown ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 79daa6d063SJed Brown if (lsc->scale) { 80daa6d063SJed Brown ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 81daa6d063SJed Brown } 82daa6d063SJed Brown ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 83daa6d063SJed Brown if (lsc->scale) { 84daa6d063SJed Brown ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 85daa6d063SJed Brown } 86daa6d063SJed Brown ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 87daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 88c0decd05SBarry Smith ierr = KSPCheckSolve(lsc->kspL,pc,y);CHKERRQ(ierr); 89daa6d063SJed Brown PetscFunctionReturn(0); 90daa6d063SJed Brown } 91daa6d063SJed Brown 92a06653b4SBarry Smith static PetscErrorCode PCReset_LSC(PC pc) 93daa6d063SJed Brown { 94daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 95daa6d063SJed Brown PetscErrorCode ierr; 96daa6d063SJed Brown 97daa6d063SJed Brown PetscFunctionBegin; 986bf464f9SBarry Smith ierr = VecDestroy(&lsc->x0);CHKERRQ(ierr); 996bf464f9SBarry Smith ierr = VecDestroy(&lsc->y0);CHKERRQ(ierr); 1006bf464f9SBarry Smith ierr = VecDestroy(&lsc->x1);CHKERRQ(ierr); 1016bf464f9SBarry Smith ierr = VecDestroy(&lsc->scale);CHKERRQ(ierr); 1026bf464f9SBarry Smith ierr = KSPDestroy(&lsc->kspL);CHKERRQ(ierr); 1036bf464f9SBarry Smith ierr = MatDestroy(&lsc->L);CHKERRQ(ierr); 104a06653b4SBarry Smith PetscFunctionReturn(0); 105a06653b4SBarry Smith } 106a06653b4SBarry Smith 107a06653b4SBarry Smith static PetscErrorCode PCDestroy_LSC(PC pc) 108a06653b4SBarry Smith { 109a06653b4SBarry Smith PetscErrorCode ierr; 110a06653b4SBarry Smith 111a06653b4SBarry Smith PetscFunctionBegin; 112a06653b4SBarry Smith ierr = PCReset_LSC(pc);CHKERRQ(ierr); 113c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 114daa6d063SJed Brown PetscFunctionReturn(0); 115daa6d063SJed Brown } 116daa6d063SJed Brown 1174416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LSC(PetscOptionItems *PetscOptionsObject,PC pc) 118daa6d063SJed Brown { 119daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 120daa6d063SJed Brown PetscErrorCode ierr; 121daa6d063SJed Brown 122daa6d063SJed Brown PetscFunctionBegin; 123e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"LSC options");CHKERRQ(ierr); 124e55864a3SBarry Smith { 1250298fd71SBarry Smith ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,NULL);CHKERRQ(ierr); 1261a1499c8SBarry Smith } 127daa6d063SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 128daa6d063SJed Brown PetscFunctionReturn(0); 129daa6d063SJed Brown } 130daa6d063SJed Brown 1317e8cb189SBarry Smith static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 1327e8cb189SBarry Smith { 1337e8cb189SBarry Smith PC_LSC *jac = (PC_LSC*)pc->data; 1347e8cb189SBarry Smith PetscErrorCode ierr; 1357e8cb189SBarry Smith PetscBool iascii; 1367e8cb189SBarry Smith 1377e8cb189SBarry Smith PetscFunctionBegin; 138251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1397e8cb189SBarry Smith if (iascii) { 1407e8cb189SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 141ca0c957dSBarry Smith if (jac->kspL) { 1427e8cb189SBarry Smith ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 143ca0c957dSBarry Smith } else { 144ca0c957dSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PCLSC KSP object not yet created, hence cannot display");CHKERRQ(ierr); 145ca0c957dSBarry Smith } 1467e8cb189SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 14711aeaf0aSBarry Smith } 1487e8cb189SBarry Smith PetscFunctionReturn(0); 1497e8cb189SBarry Smith } 1507e8cb189SBarry Smith 151daa6d063SJed Brown /*MC 152daa6d063SJed Brown PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 153daa6d063SJed Brown 154daa6d063SJed Brown Options Database Key: 155daa6d063SJed Brown . -pc_lsc_scale_diag - Use the diagonal of A for scaling 156daa6d063SJed Brown 157daa6d063SJed Brown Level: intermediate 158daa6d063SJed Brown 159daa6d063SJed Brown Notes: 160daa6d063SJed Brown This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 161daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 162daa6d063SJed Brown 163daa6d063SJed Brown .vb 1647e8cb189SBarry Smith S = A11 - A10 inv(A00) A01 165daa6d063SJed Brown .ve 166daa6d063SJed Brown 1677e8cb189SBarry Smith PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 168daa6d063SJed Brown inv(S) is given by 169daa6d063SJed Brown 170daa6d063SJed Brown .vb 1717e8cb189SBarry Smith inv(A10 A01) A10 A00 A01 inv(A10 A01) 172daa6d063SJed Brown .ve 173daa6d063SJed Brown 1747e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 1751c6f1693SPatrick Sanan usually more efficient anyway). In the case of incompressible flow, A10 A01 is a Laplacian; call it L. The current 176daa6d063SJed Brown interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 177daa6d063SJed Brown 17862279869SBarry Smith If you had called KSPSetOperators(ksp,S,Sp), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 179daa6d063SJed Brown like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 180daa6d063SJed Brown For example, you might have setup code like this 181daa6d063SJed Brown 182daa6d063SJed Brown .vb 183daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 184daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 185daa6d063SJed Brown .ve 186daa6d063SJed Brown 187daa6d063SJed Brown And then your Jacobian assembly would look like 188daa6d063SJed Brown 189daa6d063SJed Brown .vb 190daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 191daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 192daa6d063SJed Brown if (L) { assembly L } 193daa6d063SJed Brown if (Lp) { assemble Lp } 194daa6d063SJed Brown .ve 195daa6d063SJed Brown 196daa6d063SJed Brown With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 197daa6d063SJed Brown 198daa6d063SJed Brown .vb 199daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 200daa6d063SJed Brown .ve 201daa6d063SJed Brown 202daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 203daa6d063SJed Brown 204a674a499SJed Brown References: 205*606c0280SSatish Balay + * - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006. 206*606c0280SSatish Balay - * - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001. 207a674a499SJed Brown 208daa6d063SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 20929f8a81cSJed Brown PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(), 210daa6d063SJed Brown MatCreateSchurComplement() 211daa6d063SJed Brown M*/ 212daa6d063SJed Brown 2138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc) 214daa6d063SJed Brown { 215daa6d063SJed Brown PC_LSC *lsc; 216daa6d063SJed Brown PetscErrorCode ierr; 217daa6d063SJed Brown 218daa6d063SJed Brown PetscFunctionBegin; 219b00a9115SJed Brown ierr = PetscNewLog(pc,&lsc);CHKERRQ(ierr); 220daa6d063SJed Brown pc->data = (void*)lsc; 221daa6d063SJed Brown 222daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 2230a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 224daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 225a06653b4SBarry Smith pc->ops->reset = PCReset_LSC; 226daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 227daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 2287e8cb189SBarry Smith pc->ops->view = PCView_LSC; 2290a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 230daa6d063SJed Brown PetscFunctionReturn(0); 231daa6d063SJed Brown } 232