1daa6d063SJed Brown 2c6db04a5SJed Brown #include <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); 23daa6d063SJed Brown ierr = KSPCreate(((PetscObject)pc)->comm,&lsc->kspL);CHKERRQ(ierr); 247e8cb189SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);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; 447e8cb189SBarry Smith Mat L,Lp,B,C; 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); 517e8cb189SBarry Smith if (!L) { 527e8cb189SBarry Smith ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr); 537e8cb189SBarry Smith if (!lsc->L) { 547e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 557e8cb189SBarry Smith } else { 567e8cb189SBarry Smith ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 577e8cb189SBarry Smith } 587e8cb189SBarry Smith Lp = L = lsc->L; 597e8cb189SBarry Smith } 60daa6d063SJed Brown if (lsc->scale) { 61daa6d063SJed Brown Mat Ap; 62daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,&Ap,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 63daa6d063SJed Brown ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 64daa6d063SJed Brown ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 65daa6d063SJed Brown } 66daa6d063SJed Brown ierr = KSPSetOperators(lsc->kspL,L,Lp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 67daa6d063SJed Brown PetscFunctionReturn(0); 68daa6d063SJed Brown } 69daa6d063SJed Brown 70daa6d063SJed Brown #undef __FUNCT__ 71daa6d063SJed Brown #define __FUNCT__ "PCApply_LSC" 72daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 73daa6d063SJed Brown { 74daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 75daa6d063SJed Brown Mat A,B,C; 76daa6d063SJed Brown PetscErrorCode ierr; 77daa6d063SJed Brown 78daa6d063SJed Brown PetscFunctionBegin; 79daa6d063SJed Brown ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr); 80daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 81daa6d063SJed Brown ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 82daa6d063SJed Brown if (lsc->scale) { 83daa6d063SJed Brown ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 84daa6d063SJed Brown } 85daa6d063SJed Brown ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 86daa6d063SJed Brown if (lsc->scale) { 87daa6d063SJed Brown ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 88daa6d063SJed Brown } 89daa6d063SJed Brown ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 90daa6d063SJed Brown ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 91daa6d063SJed Brown PetscFunctionReturn(0); 92daa6d063SJed Brown } 93daa6d063SJed Brown 94daa6d063SJed Brown #undef __FUNCT__ 95a06653b4SBarry Smith #define __FUNCT__ "PCReset_LSC" 96a06653b4SBarry Smith static PetscErrorCode PCReset_LSC(PC pc) 97daa6d063SJed Brown { 98daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 99daa6d063SJed Brown PetscErrorCode ierr; 100daa6d063SJed Brown 101daa6d063SJed Brown PetscFunctionBegin; 102*6bf464f9SBarry Smith ierr = VecDestroy(&lsc->x0);CHKERRQ(ierr); 103*6bf464f9SBarry Smith ierr = VecDestroy(&lsc->y0);CHKERRQ(ierr); 104*6bf464f9SBarry Smith ierr = VecDestroy(&lsc->x1);CHKERRQ(ierr); 105*6bf464f9SBarry Smith ierr = VecDestroy(&lsc->scale);CHKERRQ(ierr); 106*6bf464f9SBarry Smith ierr = KSPDestroy(&lsc->kspL);CHKERRQ(ierr); 107*6bf464f9SBarry Smith ierr = MatDestroy(&lsc->L);CHKERRQ(ierr); 108a06653b4SBarry Smith PetscFunctionReturn(0); 109a06653b4SBarry Smith } 110a06653b4SBarry Smith 111a06653b4SBarry Smith #undef __FUNCT__ 112a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_LSC" 113a06653b4SBarry Smith static PetscErrorCode PCDestroy_LSC(PC pc) 114a06653b4SBarry Smith { 115a06653b4SBarry Smith PetscErrorCode ierr; 116a06653b4SBarry Smith 117a06653b4SBarry Smith PetscFunctionBegin; 118a06653b4SBarry Smith ierr = PCReset_LSC(pc);CHKERRQ(ierr); 119c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 120daa6d063SJed Brown PetscFunctionReturn(0); 121daa6d063SJed Brown } 122daa6d063SJed Brown 123daa6d063SJed Brown #undef __FUNCT__ 124daa6d063SJed Brown #define __FUNCT__ "PCSetFromOptions_LSC" 125daa6d063SJed Brown static PetscErrorCode PCSetFromOptions_LSC(PC pc) 126daa6d063SJed Brown { 127daa6d063SJed Brown PC_LSC *lsc = (PC_LSC*)pc->data; 128daa6d063SJed Brown PetscErrorCode ierr; 129daa6d063SJed Brown 130daa6d063SJed Brown PetscFunctionBegin; 131daa6d063SJed Brown ierr = PetscOptionsHead("LSC options");CHKERRQ(ierr); 132daa6d063SJed Brown { 133acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,PETSC_NULL);CHKERRQ(ierr); 134daa6d063SJed Brown } 135daa6d063SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 136daa6d063SJed Brown PetscFunctionReturn(0); 137daa6d063SJed Brown } 138daa6d063SJed Brown 1397e8cb189SBarry Smith #undef __FUNCT__ 1407e8cb189SBarry Smith #define __FUNCT__ "PCView_LSC" 1417e8cb189SBarry Smith static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 1427e8cb189SBarry Smith { 1437e8cb189SBarry Smith PC_LSC *jac = (PC_LSC*)pc->data; 1447e8cb189SBarry Smith PetscErrorCode ierr; 1457e8cb189SBarry Smith PetscBool iascii; 1467e8cb189SBarry Smith 1477e8cb189SBarry Smith PetscFunctionBegin; 1487e8cb189SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1497e8cb189SBarry Smith if (iascii) { 1507e8cb189SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1517e8cb189SBarry Smith ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 1527e8cb189SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1537e8cb189SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for LSC",((PetscObject)viewer)->type_name); 1547e8cb189SBarry Smith PetscFunctionReturn(0); 1557e8cb189SBarry Smith } 1567e8cb189SBarry Smith 157daa6d063SJed Brown /*MC 158daa6d063SJed Brown PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 159daa6d063SJed Brown 160daa6d063SJed Brown Options Database Key: 161daa6d063SJed Brown . -pc_lsc_scale_diag - Use the diagonal of A for scaling 162daa6d063SJed Brown 163daa6d063SJed Brown Level: intermediate 164daa6d063SJed Brown 165daa6d063SJed Brown Notes: 166daa6d063SJed Brown This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 167daa6d063SJed Brown it can be used for any Schur complement system. Consider the Schur complement 168daa6d063SJed Brown 169daa6d063SJed Brown .vb 1707e8cb189SBarry Smith S = A11 - A10 inv(A00) A01 171daa6d063SJed Brown .ve 172daa6d063SJed Brown 1737e8cb189SBarry Smith PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 174daa6d063SJed Brown inv(S) is given by 175daa6d063SJed Brown 176daa6d063SJed Brown .vb 1777e8cb189SBarry Smith inv(A10 A01) A10 A00 A01 inv(A10 A01) 178daa6d063SJed Brown .ve 179daa6d063SJed Brown 1807e8cb189SBarry Smith The product A10 A01 can be computed for you, but you can provide it (this is 1817e8cb189SBarry Smith usually more efficient anyway). In the case of incompressible flow, A10 A10 is a Laplacian, call it L. The current 182daa6d063SJed Brown interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 183daa6d063SJed Brown 184daa6d063SJed Brown If you had called KSPSetOperators(ksp,S,Sp,flg), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 185daa6d063SJed Brown like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 186daa6d063SJed Brown For example, you might have setup code like this 187daa6d063SJed Brown 188daa6d063SJed Brown .vb 189daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 190daa6d063SJed Brown PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 191daa6d063SJed Brown .ve 192daa6d063SJed Brown 193daa6d063SJed Brown And then your Jacobian assembly would look like 194daa6d063SJed Brown 195daa6d063SJed Brown .vb 196daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 197daa6d063SJed Brown PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 198daa6d063SJed Brown if (L) { assembly L } 199daa6d063SJed Brown if (Lp) { assemble Lp } 200daa6d063SJed Brown .ve 201daa6d063SJed Brown 202daa6d063SJed Brown With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 203daa6d063SJed Brown 204daa6d063SJed Brown .vb 205daa6d063SJed Brown -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 206daa6d063SJed Brown .ve 207daa6d063SJed Brown 208daa6d063SJed Brown Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 209daa6d063SJed Brown 210daa6d063SJed Brown Concepts: physics based preconditioners, block preconditioners 211daa6d063SJed Brown 212daa6d063SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 213daa6d063SJed Brown PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition(), 214daa6d063SJed Brown MatCreateSchurComplement() 215daa6d063SJed Brown M*/ 216daa6d063SJed Brown 217daa6d063SJed Brown EXTERN_C_BEGIN 218daa6d063SJed Brown #undef __FUNCT__ 219daa6d063SJed Brown #define __FUNCT__ "PCCreate_LSC" 2207087cfbeSBarry Smith PetscErrorCode PCCreate_LSC(PC pc) 221daa6d063SJed Brown { 222daa6d063SJed Brown PC_LSC *lsc; 223daa6d063SJed Brown PetscErrorCode ierr; 224daa6d063SJed Brown 225daa6d063SJed Brown PetscFunctionBegin; 226daa6d063SJed Brown ierr = PetscNewLog(pc,PC_LSC,&lsc);CHKERRQ(ierr); 227daa6d063SJed Brown pc->data = (void*)lsc; 228daa6d063SJed Brown 229daa6d063SJed Brown pc->ops->apply = PCApply_LSC; 230daa6d063SJed Brown pc->ops->applytranspose = 0; 231daa6d063SJed Brown pc->ops->setup = PCSetUp_LSC; 232a06653b4SBarry Smith pc->ops->reset = PCReset_LSC; 233daa6d063SJed Brown pc->ops->destroy = PCDestroy_LSC; 234daa6d063SJed Brown pc->ops->setfromoptions = PCSetFromOptions_LSC; 2357e8cb189SBarry Smith pc->ops->view = PCView_LSC; 236daa6d063SJed Brown pc->ops->applyrichardson = 0; 237daa6d063SJed Brown PetscFunctionReturn(0); 238daa6d063SJed Brown } 239daa6d063SJed Brown EXTERN_C_END 240