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; 12daa6d063SJed Brown } PC_LSC; 13daa6d063SJed Brown 14daa6d063SJed Brown #undef __FUNCT__ 15daa6d063SJed Brown #define __FUNCT__ "PCLSCAllocate_Private" 16daa6d063SJed Brown static PetscErrorCode PCLSCAllocate_Private(PC pc) 17daa6d063SJed Brown { 18daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 19daa6d063SJed Brown Mat A; 20daa6d063SJed Brown PetscErrorCode ierr; 21daa6d063SJed Brown 22daa6d063SJed Brown PetscFunctionBegin; 23daa6d063SJed Brown if (lsc->allocated) PetscFunctionReturn(0); 24daa6d063SJed Brown ierr = KSPCreate(((PetscObject)pc)->comm,&lsc->kspL);CHKERRQ(ierr); 25daa6d063SJed Brown ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr); 26daa6d063SJed Brown ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr); 27daa6d063SJed Brown ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr); 28daa6d063SJed Brown ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr); 29daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 30daa6d063SJed Brown ierr = MatGetVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr); 31daa6d063SJed Brown ierr = MatGetVecs(pc->pmat,&lsc->x1,PETSC_NULL);CHKERRQ(ierr); 32daa6d063SJed Brown if (lsc->scalediag) { 33daa6d063SJed Brown ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr); 34daa6d063SJed Brown } 35daa6d063SJed Brown lsc->allocated = PETSC_TRUE; 36daa6d063SJed Brown PetscFunctionReturn(0); 37daa6d063SJed Brown } 38daa6d063SJed Brown 39daa6d063SJed Brown #undef __FUNCT__ 40daa6d063SJed Brown #define __FUNCT__ "PCSetUp_LSC" 41daa6d063SJed Brown static PetscErrorCode PCSetUp_LSC(PC pc) 42daa6d063SJed Brown { 43daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 44daa6d063SJed Brown Mat L,Lp; 45daa6d063SJed Brown PetscErrorCode ierr; 46daa6d063SJed Brown 47daa6d063SJed Brown PetscFunctionBegin; 48daa6d063SJed Brown ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr); 49daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr); 50daa6d063SJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr); 51daa6d063SJed Brown if (lsc->scale) { 52daa6d063SJed Brown Mat Ap; 53daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,&Ap,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 54daa6d063SJed Brown ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 55daa6d063SJed Brown ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 56daa6d063SJed Brown } 57daa6d063SJed Brown ierr = KSPSetOperators(lsc->kspL,L,Lp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 58daa6d063SJed Brown PetscFunctionReturn(0); 59daa6d063SJed Brown } 60daa6d063SJed Brown 61daa6d063SJed Brown #undef __FUNCT__ 62daa6d063SJed Brown #define __FUNCT__ "PCApply_LSC" 63daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 64daa6d063SJed Brown { 65daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 66daa6d063SJed Brown Mat A,B,C; 67daa6d063SJed Brown PetscErrorCode ierr; 68daa6d063SJed Brown 69daa6d063SJed Brown PetscFunctionBegin; 70daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr); 71daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 72daa6d063SJed Brown ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 73daa6d063SJed Brown if (lsc->scale) { 74daa6d063SJed Brown ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 75daa6d063SJed Brown } 76daa6d063SJed Brown ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 77daa6d063SJed Brown if (lsc->scale) { 78daa6d063SJed Brown ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 79daa6d063SJed Brown } 80daa6d063SJed Brown ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 81daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 82daa6d063SJed Brown PetscFunctionReturn(0); 83daa6d063SJed Brown } 84daa6d063SJed Brown 85daa6d063SJed Brown #undef __FUNCT__ 86daa6d063SJed Brown #define __FUNCT__ "PCDestroy_LSC" 87daa6d063SJed Brown static PetscErrorCode PCDestroy_LSC(PC pc) 88daa6d063SJed Brown { 89daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 90daa6d063SJed Brown PetscErrorCode ierr; 91daa6d063SJed Brown 92daa6d063SJed Brown PetscFunctionBegin; 93daa6d063SJed Brown if (lsc->x0) {ierr = VecDestroy(lsc->x0);CHKERRQ(ierr);} 94daa6d063SJed Brown if (lsc->y0) {ierr = VecDestroy(lsc->y0);CHKERRQ(ierr);} 95daa6d063SJed Brown if (lsc->x1) {ierr = VecDestroy(lsc->x1);CHKERRQ(ierr);} 96daa6d063SJed Brown if (lsc->scale) {ierr = VecDestroy(lsc->scale);CHKERRQ(ierr);} 97daa6d063SJed Brown if (lsc->kspL) {ierr = KSPDestroy(lsc->kspL);CHKERRQ(ierr);} 98daa6d063SJed Brown ierr = PetscFree(lsc);CHKERRQ(ierr); 99daa6d063SJed Brown PetscFunctionReturn(0); 100daa6d063SJed Brown } 101daa6d063SJed Brown 102daa6d063SJed Brown #undef __FUNCT__ 103daa6d063SJed Brown #define __FUNCT__ "PCSetFromOptions_LSC" 104daa6d063SJed Brown static PetscErrorCode PCSetFromOptions_LSC(PC pc) 105daa6d063SJed Brown { 106daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 107daa6d063SJed Brown PetscErrorCode ierr; 108daa6d063SJed Brown 109daa6d063SJed Brown PetscFunctionBegin; 110daa6d063SJed Brown ierr = PetscOptionsHead("LSC options");CHKERRQ(ierr); 111daa6d063SJed Brown { 112*acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,PETSC_NULL);CHKERRQ(ierr); 113daa6d063SJed Brown } 114daa6d063SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 115daa6d063SJed Brown PetscFunctionReturn(0); 116daa6d063SJed Brown } 117daa6d063SJed Brown 118daa6d063SJed Brown /*MC 119daa6d063SJed Brown PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 120daa6d063SJed Brown 121daa6d063SJed Brown Options Database Key: 122daa6d063SJed Brown . -pc_lsc_scale_diag - Use the diagonal of A for scaling 123daa6d063SJed Brown 124daa6d063SJed Brown Level: intermediate 125daa6d063SJed Brown 126daa6d063SJed Brown Notes: 127daa6d063SJed Brown This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 128daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 129daa6d063SJed Brown 130daa6d063SJed Brown .vb 131daa6d063SJed Brown S = D - C inv(A) B 132daa6d063SJed Brown .ve 133daa6d063SJed Brown 134daa6d063SJed Brown PCLSC currently doesn't do anything with D, so let's assume it is 0. The idea is that a good approximation to 135daa6d063SJed Brown inv(S) is given by 136daa6d063SJed Brown 137daa6d063SJed Brown .vb 138daa6d063SJed Brown inv(CB) C A B inv(CB) 139daa6d063SJed Brown .ve 140daa6d063SJed Brown 141daa6d063SJed Brown At some point, we'll be able to form the product CB for you, but for now the application has to provide it (this is 142daa6d063SJed Brown usually more efficient anyway). In the case of incompressible flow, CB is a Laplacian, call it L. The current 143daa6d063SJed Brown interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 144daa6d063SJed Brown 145daa6d063SJed Brown If you had called KSPSetOperators(ksp,S,Sp,flg), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 146daa6d063SJed Brown like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 147daa6d063SJed Brown For example, you might have setup code like this 148daa6d063SJed Brown 149daa6d063SJed Brown .vb 150daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 151daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 152daa6d063SJed Brown .ve 153daa6d063SJed Brown 154daa6d063SJed Brown And then your Jacobian assembly would look like 155daa6d063SJed Brown 156daa6d063SJed Brown .vb 157daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 158daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 159daa6d063SJed Brown if (L) { assembly L } 160daa6d063SJed Brown if (Lp) { assemble Lp } 161daa6d063SJed Brown .ve 162daa6d063SJed Brown 163daa6d063SJed Brown With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 164daa6d063SJed Brown 165daa6d063SJed Brown .vb 166daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 167daa6d063SJed Brown .ve 168daa6d063SJed Brown 169daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 170daa6d063SJed Brown 171daa6d063SJed Brown Concepts: physics based preconditioners, block preconditioners 172daa6d063SJed Brown 173daa6d063SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 174daa6d063SJed Brown PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition(), 175daa6d063SJed Brown MatCreateSchurComplement() 176daa6d063SJed Brown M*/ 177daa6d063SJed Brown 178daa6d063SJed Brown EXTERN_C_BEGIN 179daa6d063SJed Brown #undef __FUNCT__ 180daa6d063SJed Brown #define __FUNCT__ "PCCreate_LSC" 181daa6d063SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LSC(PC pc) 182daa6d063SJed Brown { 183daa6d063SJed Brown PC_LSC *lsc; 184daa6d063SJed Brown PetscErrorCode ierr; 185daa6d063SJed Brown 186daa6d063SJed Brown PetscFunctionBegin; 187daa6d063SJed Brown ierr = PetscNewLog(pc,PC_LSC,&lsc);CHKERRQ(ierr); 188daa6d063SJed Brown pc->data = (void*)lsc; 189daa6d063SJed Brown 190daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 191daa6d063SJed Brown pc->ops->applytranspose = 0; 192daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 193daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 194daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 195daa6d063SJed Brown pc->ops->view = 0; 196daa6d063SJed Brown pc->ops->applyrichardson = 0; 197daa6d063SJed Brown PetscFunctionReturn(0); 198daa6d063SJed Brown } 199daa6d063SJed Brown EXTERN_C_END 200