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