1 2 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 3 4 typedef struct { 5 PetscBool allocated; 6 PetscBool scalediag; 7 KSP kspL; 8 Vec scale; 9 Vec x0,y0,x1; 10 Mat L; /* keep a copy to reuse when obtained with L = A10*A01 */ 11 } PC_LSC; 12 13 static PetscErrorCode PCLSCAllocate_Private(PC pc) 14 { 15 PC_LSC *lsc = (PC_LSC*)pc->data; 16 Mat A; 17 PetscErrorCode ierr; 18 19 PetscFunctionBegin; 20 if (lsc->allocated) PetscFunctionReturn(0); 21 ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&lsc->kspL);CHKERRQ(ierr); 22 ierr = KSPSetErrorIfNotConverged(lsc->kspL,pc->erroriffailure);CHKERRQ(ierr); 23 ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);CHKERRQ(ierr); 24 ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr); 25 ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr); 26 ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr); 27 ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 28 ierr = MatCreateVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr); 29 ierr = MatCreateVecs(pc->pmat,&lsc->x1,NULL);CHKERRQ(ierr); 30 if (lsc->scalediag) { 31 ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr); 32 } 33 lsc->allocated = PETSC_TRUE; 34 PetscFunctionReturn(0); 35 } 36 37 static PetscErrorCode PCSetUp_LSC(PC pc) 38 { 39 PC_LSC *lsc = (PC_LSC*)pc->data; 40 Mat L,Lp,B,C; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr); 45 ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr); 46 if (!L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr);} 47 ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr); 48 if (!Lp) {ierr = PetscObjectQuery((PetscObject)pc->mat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr);} 49 if (!L) { 50 ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr); 51 if (!lsc->L) { 52 ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 53 } else { 54 ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr); 55 } 56 Lp = L = lsc->L; 57 } 58 if (lsc->scale) { 59 Mat Ap; 60 ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,&Ap,NULL,NULL,NULL);CHKERRQ(ierr); 61 ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */ 62 ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr); 63 } 64 ierr = KSPSetOperators(lsc->kspL,L,Lp);CHKERRQ(ierr); 65 ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr); 66 PetscFunctionReturn(0); 67 } 68 69 static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y) 70 { 71 PC_LSC *lsc = (PC_LSC*)pc->data; 72 Mat A,B,C; 73 PetscErrorCode ierr; 74 75 PetscFunctionBegin; 76 ierr = MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,&B,&C,NULL);CHKERRQ(ierr); 77 ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr); 78 ierr = KSPCheckSolve(lsc->kspL,pc,lsc->x1);CHKERRQ(ierr); 79 ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 80 if (lsc->scale) { 81 ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 82 } 83 ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 84 if (lsc->scale) { 85 ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 86 } 87 ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 88 ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 89 ierr = KSPCheckSolve(lsc->kspL,pc,y);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 static PetscErrorCode PCReset_LSC(PC pc) 94 { 95 PC_LSC *lsc = (PC_LSC*)pc->data; 96 PetscErrorCode ierr; 97 98 PetscFunctionBegin; 99 ierr = VecDestroy(&lsc->x0);CHKERRQ(ierr); 100 ierr = VecDestroy(&lsc->y0);CHKERRQ(ierr); 101 ierr = VecDestroy(&lsc->x1);CHKERRQ(ierr); 102 ierr = VecDestroy(&lsc->scale);CHKERRQ(ierr); 103 ierr = KSPDestroy(&lsc->kspL);CHKERRQ(ierr); 104 ierr = MatDestroy(&lsc->L);CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 static PetscErrorCode PCDestroy_LSC(PC pc) 109 { 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 ierr = PCReset_LSC(pc);CHKERRQ(ierr); 114 ierr = PetscFree(pc->data);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 static PetscErrorCode PCSetFromOptions_LSC(PetscOptionItems *PetscOptionsObject,PC pc) 119 { 120 PC_LSC *lsc = (PC_LSC*)pc->data; 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 ierr = PetscOptionsHead(PetscOptionsObject,"LSC options");CHKERRQ(ierr); 125 { 126 ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,NULL);CHKERRQ(ierr); 127 } 128 ierr = PetscOptionsTail();CHKERRQ(ierr); 129 PetscFunctionReturn(0); 130 } 131 132 static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 133 { 134 PC_LSC *jac = (PC_LSC*)pc->data; 135 PetscErrorCode ierr; 136 PetscBool iascii; 137 138 PetscFunctionBegin; 139 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 140 if (iascii) { 141 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 142 ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 143 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 144 } 145 PetscFunctionReturn(0); 146 } 147 148 /*MC 149 PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 150 151 Options Database Key: 152 . -pc_lsc_scale_diag - Use the diagonal of A for scaling 153 154 Level: intermediate 155 156 Notes: 157 This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 158 it can be used for any Schur complement system. Consider the Schur complement 159 160 .vb 161 S = A11 - A10 inv(A00) A01 162 .ve 163 164 PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 165 inv(S) is given by 166 167 .vb 168 inv(A10 A01) A10 A00 A01 inv(A10 A01) 169 .ve 170 171 The product A10 A01 can be computed for you, but you can provide it (this is 172 usually more efficient anyway). In the case of incompressible flow, A10 A10 is a Laplacian, call it L. The current 173 interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 174 175 If you had called KSPSetOperators(ksp,S,Sp), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 176 like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 177 For example, you might have setup code like this 178 179 .vb 180 PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 181 PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 182 .ve 183 184 And then your Jacobian assembly would look like 185 186 .vb 187 PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 188 PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 189 if (L) { assembly L } 190 if (Lp) { assemble Lp } 191 .ve 192 193 With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 194 195 .vb 196 -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 197 .ve 198 199 Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 200 201 References: 202 + 1. - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006. 203 - 2. - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001. 204 205 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT, 206 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(), 207 MatCreateSchurComplement() 208 M*/ 209 210 PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc) 211 { 212 PC_LSC *lsc; 213 PetscErrorCode ierr; 214 215 PetscFunctionBegin; 216 ierr = PetscNewLog(pc,&lsc);CHKERRQ(ierr); 217 pc->data = (void*)lsc; 218 219 pc->ops->apply = PCApply_LSC; 220 pc->ops->applytranspose = 0; 221 pc->ops->setup = PCSetUp_LSC; 222 pc->ops->reset = PCReset_LSC; 223 pc->ops->destroy = PCDestroy_LSC; 224 pc->ops->setfromoptions = PCSetFromOptions_LSC; 225 pc->ops->view = PCView_LSC; 226 pc->ops->applyrichardson = 0; 227 PetscFunctionReturn(0); 228 } 229