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 = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr); 79 if (lsc->scale) { 80 ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr); 81 } 82 ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr); 83 if (lsc->scale) { 84 ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr); 85 } 86 ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr); 87 ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 91 static PetscErrorCode PCReset_LSC(PC pc) 92 { 93 PC_LSC *lsc = (PC_LSC*)pc->data; 94 PetscErrorCode ierr; 95 96 PetscFunctionBegin; 97 ierr = VecDestroy(&lsc->x0);CHKERRQ(ierr); 98 ierr = VecDestroy(&lsc->y0);CHKERRQ(ierr); 99 ierr = VecDestroy(&lsc->x1);CHKERRQ(ierr); 100 ierr = VecDestroy(&lsc->scale);CHKERRQ(ierr); 101 ierr = KSPDestroy(&lsc->kspL);CHKERRQ(ierr); 102 ierr = MatDestroy(&lsc->L);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 106 static PetscErrorCode PCDestroy_LSC(PC pc) 107 { 108 PetscErrorCode ierr; 109 110 PetscFunctionBegin; 111 ierr = PCReset_LSC(pc);CHKERRQ(ierr); 112 ierr = PetscFree(pc->data);CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 static PetscErrorCode PCSetFromOptions_LSC(PetscOptionItems *PetscOptionsObject,PC pc) 117 { 118 PC_LSC *lsc = (PC_LSC*)pc->data; 119 PetscErrorCode ierr; 120 121 PetscFunctionBegin; 122 ierr = PetscOptionsHead(PetscOptionsObject,"LSC options");CHKERRQ(ierr); 123 { 124 ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,NULL);CHKERRQ(ierr); 125 } 126 ierr = PetscOptionsTail();CHKERRQ(ierr); 127 PetscFunctionReturn(0); 128 } 129 130 static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer) 131 { 132 PC_LSC *jac = (PC_LSC*)pc->data; 133 PetscErrorCode ierr; 134 PetscBool iascii; 135 136 PetscFunctionBegin; 137 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 138 if (iascii) { 139 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 140 ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr); 141 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 142 } 143 PetscFunctionReturn(0); 144 } 145 146 /*MC 147 PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators 148 149 Options Database Key: 150 . -pc_lsc_scale_diag - Use the diagonal of A for scaling 151 152 Level: intermediate 153 154 Notes: 155 This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but 156 it can be used for any Schur complement system. Consider the Schur complement 157 158 .vb 159 S = A11 - A10 inv(A00) A01 160 .ve 161 162 PCLSC currently doesn't do anything with A11, so let's assume it is 0. The idea is that a good approximation to 163 inv(S) is given by 164 165 .vb 166 inv(A10 A01) A10 A00 A01 inv(A10 A01) 167 .ve 168 169 The product A10 A01 can be computed for you, but you can provide it (this is 170 usually more efficient anyway). In the case of incompressible flow, A10 A10 is a Laplacian, call it L. The current 171 interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix. 172 173 If you had called KSPSetOperators(ksp,S,Sp), S should have type MATSCHURCOMPLEMENT and Sp can be any type you 174 like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp". 175 For example, you might have setup code like this 176 177 .vb 178 PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L); 179 PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp); 180 .ve 181 182 And then your Jacobian assembly would look like 183 184 .vb 185 PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L); 186 PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp); 187 if (L) { assembly L } 188 if (Lp) { assemble Lp } 189 .ve 190 191 With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L 192 193 .vb 194 -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml 195 .ve 196 197 Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners. 198 199 References: 200 + 1. - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006. 201 - 2. - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001. 202 203 Concepts: physics based preconditioners, block preconditioners 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