1daa6d063SJed Brown #define PETSCKSP_DLL 2daa6d063SJed Brown 3daa6d063SJed Brown 4daa6d063SJed Brown #include "private/pcimpl.h" /*I "petscpc.h" I*/ 5daa6d063SJed Brown 6daa6d063SJed Brown typedef struct { 7ace3abfcSBarry Smith PetscBool allocated; 8ace3abfcSBarry Smith PetscBool scalediag; 9daa6d063SJed Brown KSP kspL; 10daa6d063SJed Brown Vec scale; 11daa6d063SJed Brown Vec x0,y0,x1; 12*7e8cb189SBarry Smith Mat L; /* keep a copy to reuse when obtained with L = A10*A01 */ 13daa6d063SJed Brown } PC_LSC; 14daa6d063SJed Brown 15daa6d063SJed Brown #undef __FUNCT__ 16daa6d063SJed Brown #define __FUNCT__ "PCLSCAllocate_Private" 17daa6d063SJed Brown static PetscErrorCode PCLSCAllocate_Private(PC pc) 18daa6d063SJed Brown { 19daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 20daa6d063SJed Brown Mat A; 21daa6d063SJed Brown PetscErrorCode ierr; 22daa6d063SJed Brown 23daa6d063SJed Brown PetscFunctionBegin; 24daa6d063SJed Brown if (lsc->allocated) PetscFunctionReturn(0); 25daa6d063SJed Brown ierr = KSPCreate(((PetscObject)pc)->comm,&lsc->kspL);CHKERRQ(ierr); 26*7e8cb189SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);CHKERRQ(ierr); 27daa6d063SJed Brown ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr); 28daa6d063SJed Brown ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr); 29daa6d063SJed Brown ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr); 30daa6d063SJed Brown ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr); 31daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 32daa6d063SJed Brown ierr = MatGetVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr); 33daa6d063SJed Brown ierr = MatGetVecs(pc->pmat,&lsc->x1,PETSC_NULL);CHKERRQ(ierr); 34daa6d063SJed Brown if (lsc->scalediag) { 35daa6d063SJed Brown ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr); 36daa6d063SJed Brown } 37daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 38daa6d063SJed Brown PetscFunctionReturn(0); 39daa6d063SJed Brown } 40daa6d063SJed Brown 41daa6d063SJed Brown #undef __FUNCT__ 42daa6d063SJed Brown #define __FUNCT__ "PCSetUp_LSC" 43daa6d063SJed Brown static PetscErrorCode PCSetUp_LSC(PC pc) 44daa6d063SJed Brown { 45daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 46*7e8cb189SBarry Smith Mat L,Lp,B,C; 47daa6d063SJed Brown PetscErrorCode ierr; 48daa6d063SJed Brown 49daa6d063SJed Brown PetscFunctionBegin; 50daa6d063SJed Brown ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr); 51daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr); 52daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr); 53*7e8cb189SBarry Smith if (!L) { 54*7e8cb189SBarry Smith ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr); 55*7e8cb189SBarry Smith if (!lsc->L) { 56*7e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 57*7e8cb189SBarry Smith } else { 58*7e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 59*7e8cb189SBarry Smith } 60*7e8cb189SBarry Smith Lp = L = lsc->L; 61*7e8cb189SBarry Smith } 62daa6d063SJed Brown if (lsc->scale) { 63daa6d063SJed Brown Mat Ap; 64daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,&Ap,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 65daa6d063SJed Brown ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 66daa6d063SJed Brown ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 67daa6d063SJed Brown } 68daa6d063SJed Brown ierr = KSPSetOperators(lsc->kspL,L,Lp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 69daa6d063SJed Brown PetscFunctionReturn(0); 70daa6d063SJed Brown } 71daa6d063SJed Brown 72daa6d063SJed Brown #undef __FUNCT__ 73daa6d063SJed Brown #define __FUNCT__ "PCApply_LSC" 74daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 75daa6d063SJed Brown { 76daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 77daa6d063SJed Brown Mat A,B,C; 78daa6d063SJed Brown PetscErrorCode ierr; 79daa6d063SJed Brown 80daa6d063SJed Brown PetscFunctionBegin; 81daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr); 82daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 83daa6d063SJed Brown ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 84daa6d063SJed Brown if (lsc->scale) { 85daa6d063SJed Brown ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 86daa6d063SJed Brown } 87daa6d063SJed Brown ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 88daa6d063SJed Brown if (lsc->scale) { 89daa6d063SJed Brown ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 90daa6d063SJed Brown } 91daa6d063SJed Brown ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 92daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 93daa6d063SJed Brown PetscFunctionReturn(0); 94daa6d063SJed Brown } 95daa6d063SJed Brown 96daa6d063SJed Brown #undef __FUNCT__ 97daa6d063SJed Brown #define __FUNCT__ "PCDestroy_LSC" 98daa6d063SJed Brown static PetscErrorCode PCDestroy_LSC(PC pc) 99daa6d063SJed Brown { 100daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 101daa6d063SJed Brown PetscErrorCode ierr; 102daa6d063SJed Brown 103daa6d063SJed Brown PetscFunctionBegin; 104daa6d063SJed Brown if (lsc->x0) {ierr = VecDestroy(lsc->x0);CHKERRQ(ierr);} 105daa6d063SJed Brown if (lsc->y0) {ierr = VecDestroy(lsc->y0);CHKERRQ(ierr);} 106daa6d063SJed Brown if (lsc->x1) {ierr = VecDestroy(lsc->x1);CHKERRQ(ierr);} 107daa6d063SJed Brown if (lsc->scale) {ierr = VecDestroy(lsc->scale);CHKERRQ(ierr);} 108daa6d063SJed Brown if (lsc->kspL) {ierr = KSPDestroy(lsc->kspL);CHKERRQ(ierr);} 109*7e8cb189SBarry Smith if (lsc->L) {ierr = MatDestroy(lsc->L);CHKERRQ(ierr);} 110daa6d063SJed Brown ierr = PetscFree(lsc);CHKERRQ(ierr); 111daa6d063SJed Brown PetscFunctionReturn(0); 112daa6d063SJed Brown } 113daa6d063SJed Brown 114daa6d063SJed Brown #undef __FUNCT__ 115daa6d063SJed Brown #define __FUNCT__ "PCSetFromOptions_LSC" 116daa6d063SJed Brown static PetscErrorCode PCSetFromOptions_LSC(PC pc) 117daa6d063SJed Brown { 118daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 119daa6d063SJed Brown PetscErrorCode ierr; 120daa6d063SJed Brown 121daa6d063SJed Brown PetscFunctionBegin; 122daa6d063SJed Brown ierr = PetscOptionsHead("LSC options");CHKERRQ(ierr); 123daa6d063SJed Brown { 124acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,PETSC_NULL);CHKERRQ(ierr); 125daa6d063SJed Brown } 126daa6d063SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 127daa6d063SJed Brown PetscFunctionReturn(0); 128daa6d063SJed Brown } 129daa6d063SJed Brown 130*7e8cb189SBarry Smith #undef __FUNCT__ 131*7e8cb189SBarry Smith #define __FUNCT__ "PCView_LSC" 132*7e8cb189SBarry Smith static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 133*7e8cb189SBarry Smith { 134*7e8cb189SBarry Smith PC_LSC *jac = (PC_LSC*)pc->data; 135*7e8cb189SBarry Smith PetscErrorCode ierr; 136*7e8cb189SBarry Smith PetscBool iascii; 137*7e8cb189SBarry Smith 138*7e8cb189SBarry Smith PetscFunctionBegin; 139*7e8cb189SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 140*7e8cb189SBarry Smith if (iascii) { 141*7e8cb189SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 142*7e8cb189SBarry Smith ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 143*7e8cb189SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 144*7e8cb189SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for LSC",((PetscObject)viewer)->type_name); 145*7e8cb189SBarry Smith PetscFunctionReturn(0); 146*7e8cb189SBarry Smith } 147*7e8cb189SBarry 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 161*7e8cb189SBarry Smith S = A11 - A10 inv(A00) A01 162daa6d063SJed Brown .ve 163daa6d063SJed Brown 164*7e8cb189SBarry 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 168*7e8cb189SBarry Smith inv(A10 A01) A10 A00 A01 inv(A10 A01) 169daa6d063SJed Brown .ve 170daa6d063SJed Brown 171*7e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 172*7e8cb189SBarry 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 175daa6d063SJed Brown If you had called KSPSetOperators(ksp,S,Sp,flg), 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 201daa6d063SJed Brown Concepts: physics based preconditioners, block preconditioners 202daa6d063SJed Brown 203daa6d063SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 204daa6d063SJed Brown PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition(), 205daa6d063SJed Brown MatCreateSchurComplement() 206daa6d063SJed Brown M*/ 207daa6d063SJed Brown 208daa6d063SJed Brown EXTERN_C_BEGIN 209daa6d063SJed Brown #undef __FUNCT__ 210daa6d063SJed Brown #define __FUNCT__ "PCCreate_LSC" 2117087cfbeSBarry Smith PetscErrorCode PCCreate_LSC(PC pc) 212daa6d063SJed Brown { 213daa6d063SJed Brown PC_LSC *lsc; 214daa6d063SJed Brown PetscErrorCode ierr; 215daa6d063SJed Brown 216daa6d063SJed Brown PetscFunctionBegin; 217daa6d063SJed Brown ierr = PetscNewLog(pc,PC_LSC,&lsc);CHKERRQ(ierr); 218daa6d063SJed Brown pc->data = (void*)lsc; 219daa6d063SJed Brown 220daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 221daa6d063SJed Brown pc->ops->applytranspose = 0; 222daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 223daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 224daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 225*7e8cb189SBarry Smith pc->ops->view = PCView_LSC; 226daa6d063SJed Brown pc->ops->applyrichardson = 0; 227daa6d063SJed Brown PetscFunctionReturn(0); 228daa6d063SJed Brown } 229daa6d063SJed Brown EXTERN_C_END 230