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