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