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 #undef __FUNCT__ 14daa6d063SJed Brown #define __FUNCT__ "PCLSCAllocate_Private" 15daa6d063SJed Brown static PetscErrorCode PCLSCAllocate_Private(PC pc) 16daa6d063SJed Brown { 17daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 18daa6d063SJed Brown Mat A; 19daa6d063SJed Brown PetscErrorCode ierr; 20daa6d063SJed Brown 21daa6d063SJed Brown PetscFunctionBegin; 22daa6d063SJed Brown if (lsc->allocated) PetscFunctionReturn(0); 23ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&lsc->kspL);CHKERRQ(ierr); 24*422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(lsc->kspL,pc->erroriffailure);CHKERRQ(ierr); 257e8cb189SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);CHKERRQ(ierr); 26daa6d063SJed Brown ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr); 27daa6d063SJed Brown ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr); 28daa6d063SJed Brown ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr); 29daa6d063SJed Brown ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr); 30bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 312a7a6963SBarry Smith ierr = MatCreateVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr); 322a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&lsc->x1,NULL);CHKERRQ(ierr); 33daa6d063SJed Brown if (lsc->scalediag) { 34daa6d063SJed Brown ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr); 35daa6d063SJed Brown } 36daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 37daa6d063SJed Brown PetscFunctionReturn(0); 38daa6d063SJed Brown } 39daa6d063SJed Brown 40daa6d063SJed Brown #undef __FUNCT__ 41daa6d063SJed Brown #define __FUNCT__ "PCSetUp_LSC" 42daa6d063SJed Brown static PetscErrorCode PCSetUp_LSC(PC pc) 43daa6d063SJed Brown { 44daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 457e8cb189SBarry Smith Mat L,Lp,B,C; 46daa6d063SJed Brown PetscErrorCode ierr; 47daa6d063SJed Brown 48daa6d063SJed Brown PetscFunctionBegin; 49daa6d063SJed Brown ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr); 50093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr); 51093c86ffSJed Brown if (!L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr);} 52daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr); 53093c86ffSJed Brown if (!Lp) {ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr);} 547e8cb189SBarry Smith if (!L) { 55bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr); 567e8cb189SBarry Smith if (!lsc->L) { 577e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 587e8cb189SBarry Smith } else { 597e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 607e8cb189SBarry Smith } 617e8cb189SBarry Smith Lp = L = lsc->L; 627e8cb189SBarry Smith } 63daa6d063SJed Brown if (lsc->scale) { 64daa6d063SJed Brown Mat Ap; 65bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,&Ap,NULL,NULL,NULL);CHKERRQ(ierr); 66daa6d063SJed Brown ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 67daa6d063SJed Brown ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 68daa6d063SJed Brown } 6923ee1639SBarry Smith ierr = KSPSetOperators(lsc->kspL,L,Lp);CHKERRQ(ierr); 70daa6d063SJed Brown PetscFunctionReturn(0); 71daa6d063SJed Brown } 72daa6d063SJed Brown 73daa6d063SJed Brown #undef __FUNCT__ 74daa6d063SJed Brown #define __FUNCT__ "PCApply_LSC" 75daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 76daa6d063SJed Brown { 77daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 78daa6d063SJed Brown Mat A,B,C; 79daa6d063SJed Brown PetscErrorCode ierr; 80daa6d063SJed Brown 81daa6d063SJed Brown PetscFunctionBegin; 82bee83525SDmitry Karpeev ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,&B,&C,NULL);CHKERRQ(ierr); 83daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 84daa6d063SJed Brown ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 85daa6d063SJed Brown if (lsc->scale) { 86daa6d063SJed Brown ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 87daa6d063SJed Brown } 88daa6d063SJed Brown ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 89daa6d063SJed Brown if (lsc->scale) { 90daa6d063SJed Brown ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 91daa6d063SJed Brown } 92daa6d063SJed Brown ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 93daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 94daa6d063SJed Brown PetscFunctionReturn(0); 95daa6d063SJed Brown } 96daa6d063SJed Brown 97daa6d063SJed Brown #undef __FUNCT__ 98a06653b4SBarry Smith #define __FUNCT__ "PCReset_LSC" 99a06653b4SBarry Smith static PetscErrorCode PCReset_LSC(PC pc) 100daa6d063SJed Brown { 101daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 102daa6d063SJed Brown PetscErrorCode ierr; 103daa6d063SJed Brown 104daa6d063SJed Brown PetscFunctionBegin; 1056bf464f9SBarry Smith ierr = VecDestroy(&lsc->x0);CHKERRQ(ierr); 1066bf464f9SBarry Smith ierr = VecDestroy(&lsc->y0);CHKERRQ(ierr); 1076bf464f9SBarry Smith ierr = VecDestroy(&lsc->x1);CHKERRQ(ierr); 1086bf464f9SBarry Smith ierr = VecDestroy(&lsc->scale);CHKERRQ(ierr); 1096bf464f9SBarry Smith ierr = KSPDestroy(&lsc->kspL);CHKERRQ(ierr); 1106bf464f9SBarry Smith ierr = MatDestroy(&lsc->L);CHKERRQ(ierr); 111a06653b4SBarry Smith PetscFunctionReturn(0); 112a06653b4SBarry Smith } 113a06653b4SBarry Smith 114a06653b4SBarry Smith #undef __FUNCT__ 115a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_LSC" 116a06653b4SBarry Smith static PetscErrorCode PCDestroy_LSC(PC pc) 117a06653b4SBarry Smith { 118a06653b4SBarry Smith PetscErrorCode ierr; 119a06653b4SBarry Smith 120a06653b4SBarry Smith PetscFunctionBegin; 121a06653b4SBarry Smith ierr = PCReset_LSC(pc);CHKERRQ(ierr); 122c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 123daa6d063SJed Brown PetscFunctionReturn(0); 124daa6d063SJed Brown } 125daa6d063SJed Brown 126daa6d063SJed Brown #undef __FUNCT__ 127daa6d063SJed Brown #define __FUNCT__ "PCSetFromOptions_LSC" 1288c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_LSC(PetscOptions *PetscOptionsObject,PC pc) 129daa6d063SJed Brown { 130daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 131daa6d063SJed Brown PetscErrorCode ierr; 132daa6d063SJed Brown 133daa6d063SJed Brown PetscFunctionBegin; 134e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"LSC options");CHKERRQ(ierr); 135e55864a3SBarry Smith { 1360298fd71SBarry Smith ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,NULL);CHKERRQ(ierr); 1371a1499c8SBarry Smith } 138daa6d063SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 139daa6d063SJed Brown PetscFunctionReturn(0); 140daa6d063SJed Brown } 141daa6d063SJed Brown 1427e8cb189SBarry Smith #undef __FUNCT__ 1437e8cb189SBarry Smith #define __FUNCT__ "PCView_LSC" 1447e8cb189SBarry Smith static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 1457e8cb189SBarry Smith { 1467e8cb189SBarry Smith PC_LSC *jac = (PC_LSC*)pc->data; 1477e8cb189SBarry Smith PetscErrorCode ierr; 1487e8cb189SBarry Smith PetscBool iascii; 1497e8cb189SBarry Smith 1507e8cb189SBarry Smith PetscFunctionBegin; 151251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1527e8cb189SBarry Smith if (iascii) { 1537e8cb189SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1547e8cb189SBarry Smith ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 1557e8cb189SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 15611aeaf0aSBarry Smith } 1577e8cb189SBarry Smith PetscFunctionReturn(0); 1587e8cb189SBarry Smith } 1597e8cb189SBarry Smith 160daa6d063SJed Brown /*MC 161daa6d063SJed Brown PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 162daa6d063SJed Brown 163daa6d063SJed Brown Options Database Key: 164daa6d063SJed Brown . -pc_lsc_scale_diag - Use the diagonal of A for scaling 165daa6d063SJed Brown 166daa6d063SJed Brown Level: intermediate 167daa6d063SJed Brown 168daa6d063SJed Brown Notes: 169daa6d063SJed Brown This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 170daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 171daa6d063SJed Brown 172daa6d063SJed Brown .vb 1737e8cb189SBarry Smith S = A11 - A10 inv(A00) A01 174daa6d063SJed Brown .ve 175daa6d063SJed Brown 1767e8cb189SBarry Smith PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 177daa6d063SJed Brown inv(S) is given by 178daa6d063SJed Brown 179daa6d063SJed Brown .vb 1807e8cb189SBarry Smith inv(A10 A01) A10 A00 A01 inv(A10 A01) 181daa6d063SJed Brown .ve 182daa6d063SJed Brown 1837e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 1847e8cb189SBarry Smith usually more efficient anyway). In the case of incompressible flow, A10 A10 is a Laplacian, call it L. The current 185daa6d063SJed Brown interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 186daa6d063SJed Brown 18762279869SBarry Smith If you had called KSPSetOperators(ksp,S,Sp), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 188daa6d063SJed Brown like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 189daa6d063SJed Brown For example, you might have setup code like this 190daa6d063SJed Brown 191daa6d063SJed Brown .vb 192daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 193daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 194daa6d063SJed Brown .ve 195daa6d063SJed Brown 196daa6d063SJed Brown And then your Jacobian assembly would look like 197daa6d063SJed Brown 198daa6d063SJed Brown .vb 199daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 200daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 201daa6d063SJed Brown if (L) { assembly L } 202daa6d063SJed Brown if (Lp) { assemble Lp } 203daa6d063SJed Brown .ve 204daa6d063SJed Brown 205daa6d063SJed Brown With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 206daa6d063SJed Brown 207daa6d063SJed Brown .vb 208daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 209daa6d063SJed Brown .ve 210daa6d063SJed Brown 211daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 212daa6d063SJed Brown 213a674a499SJed Brown References: 214a674a499SJed Brown + Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006. 215a674a499SJed Brown - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier-Stokes equations for incompressible flow, 2001. 216a674a499SJed Brown 217daa6d063SJed Brown Concepts: physics based preconditioners, block preconditioners 218daa6d063SJed Brown 219daa6d063SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 22029f8a81cSJed Brown PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(), 221daa6d063SJed Brown MatCreateSchurComplement() 222daa6d063SJed Brown M*/ 223daa6d063SJed Brown 224daa6d063SJed Brown #undef __FUNCT__ 225daa6d063SJed Brown #define __FUNCT__ "PCCreate_LSC" 2268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc) 227daa6d063SJed Brown { 228daa6d063SJed Brown PC_LSC *lsc; 229daa6d063SJed Brown PetscErrorCode ierr; 230daa6d063SJed Brown 231daa6d063SJed Brown PetscFunctionBegin; 232b00a9115SJed Brown ierr = PetscNewLog(pc,&lsc);CHKERRQ(ierr); 233daa6d063SJed Brown pc->data = (void*)lsc; 234daa6d063SJed Brown 235daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 236daa6d063SJed Brown pc->ops->applytranspose = 0; 237daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 238a06653b4SBarry Smith pc->ops->reset = PCReset_LSC; 239daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 240daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 2417e8cb189SBarry Smith pc->ops->view = PCView_LSC; 242daa6d063SJed Brown pc->ops->applyrichardson = 0; 243daa6d063SJed Brown PetscFunctionReturn(0); 244daa6d063SJed Brown } 245