xref: /petsc/src/ksp/pc/impls/lsc/lsc.c (revision 9fd865de364bf685f9f753bc72082bc51a6e6586)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
2daa6d063SJed Brown 
3daa6d063SJed Brown typedef struct {
4*9fd865deSAlexander   PetscBool allocated, commute, scalediag;
5*9fd865deSAlexander   KSP       kspL, kspMass;
6*9fd865deSAlexander   Vec       Avec0, Avec1, Svec0, scale;
7*9fd865deSAlexander   Mat       L;
8daa6d063SJed Brown } PC_LSC;
9daa6d063SJed Brown 
10d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCLSCAllocate_Private(PC pc)
11d71ae5a4SJacob Faibussowitsch {
12daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
13daa6d063SJed Brown   Mat     A;
14daa6d063SJed Brown 
15daa6d063SJed Brown   PetscFunctionBegin;
163ba16761SJacob Faibussowitsch   if (lsc->allocated) PetscFunctionReturn(PETSC_SUCCESS);
179566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspL));
183821be0aSBarry Smith   PetscCall(KSPSetNestLevel(lsc->kspL, pc->kspnestlevel));
199566063dSJacob Faibussowitsch   PetscCall(KSPSetErrorIfNotConverged(lsc->kspL, pc->erroriffailure));
209566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lsc->kspL, (PetscObject)pc, 1));
219566063dSJacob Faibussowitsch   PetscCall(KSPSetType(lsc->kspL, KSPPREONLY));
229566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(lsc->kspL, ((PetscObject)pc)->prefix));
239566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(lsc->kspL, "lsc_"));
249566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, NULL, NULL, NULL));
25*9fd865deSAlexander   PetscCall(MatCreateVecs(A, &lsc->Avec0, &lsc->Avec1));
26*9fd865deSAlexander   PetscCall(MatCreateVecs(pc->pmat, &lsc->Svec0, NULL));
27*9fd865deSAlexander   if (lsc->scalediag) PetscCall(VecDuplicate(lsc->Avec0, &lsc->scale));
28*9fd865deSAlexander 
29*9fd865deSAlexander   if (lsc->commute) {
30*9fd865deSAlexander     PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspMass));
31*9fd865deSAlexander     PetscCall(KSPSetErrorIfNotConverged(lsc->kspMass, pc->erroriffailure));
32*9fd865deSAlexander     PetscCall(PetscObjectIncrementTabLevel((PetscObject)lsc->kspMass, (PetscObject)pc, 1));
33*9fd865deSAlexander     PetscCall(KSPSetType(lsc->kspMass, KSPPREONLY));
34*9fd865deSAlexander     PetscCall(KSPSetOptionsPrefix(lsc->kspMass, ((PetscObject)pc)->prefix));
35*9fd865deSAlexander     PetscCall(KSPAppendOptionsPrefix(lsc->kspMass, "lsc_mass_"));
36*9fd865deSAlexander   } else lsc->kspMass = NULL;
37*9fd865deSAlexander 
38daa6d063SJed Brown   lsc->allocated = PETSC_TRUE;
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40daa6d063SJed Brown }
41daa6d063SJed Brown 
42d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LSC(PC pc)
43d71ae5a4SJacob Faibussowitsch {
44daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
45*9fd865deSAlexander   Mat     L, Lp, Qscale;
46daa6d063SJed Brown 
47daa6d063SJed Brown   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(PCLSCAllocate_Private(pc));
49*9fd865deSAlexander 
50*9fd865deSAlexander   /* Query for L operator */
519566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_L", (PetscObject *)&L));
529566063dSJacob Faibussowitsch   if (!L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_L", (PetscObject *)&L));
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Lp", (PetscObject *)&Lp));
549566063dSJacob Faibussowitsch   if (!Lp) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Lp", (PetscObject *)&Lp));
55*9fd865deSAlexander 
56*9fd865deSAlexander   /* Query for mass operator */
57*9fd865deSAlexander   PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Qscale", (PetscObject *)&Qscale));
58*9fd865deSAlexander   if (!Qscale) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Qscale", (PetscObject *)&Qscale));
59*9fd865deSAlexander 
60*9fd865deSAlexander   if (lsc->commute) {
61*9fd865deSAlexander     PetscCheck(L || Lp, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "The user must provide an L operator for LSC preconditioning when commuting");
62*9fd865deSAlexander     if (!L && Lp) L = Lp;
63*9fd865deSAlexander     else if (L && !Lp) Lp = L;
64*9fd865deSAlexander 
65*9fd865deSAlexander     PetscCheck(Qscale, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "The user must provide a Qscale matrix for LSC preconditioning when commuting");
667e8cb189SBarry Smith   } else {
67*9fd865deSAlexander     if (lsc->scale) {
68*9fd865deSAlexander       if (!Qscale) PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, &Qscale, NULL, NULL, NULL));
69*9fd865deSAlexander       PetscCall(MatGetDiagonal(Qscale, lsc->scale));
70*9fd865deSAlexander       PetscCall(VecReciprocal(lsc->scale));
71*9fd865deSAlexander     }
72*9fd865deSAlexander     if (!L) {
73*9fd865deSAlexander       Mat B, C;
74*9fd865deSAlexander       PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, NULL, &B, &C, NULL));
75*9fd865deSAlexander       if (lsc->scale) {
76*9fd865deSAlexander         Mat CAdiaginv;
77*9fd865deSAlexander         PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &CAdiaginv));
78*9fd865deSAlexander         PetscCall(MatDiagonalScale(CAdiaginv, NULL, lsc->scale));
79*9fd865deSAlexander         if (!lsc->L) PetscCall(MatMatMult(CAdiaginv, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &lsc->L));
80*9fd865deSAlexander         else PetscCall(MatMatMult(CAdiaginv, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &lsc->L));
81*9fd865deSAlexander         PetscCall(MatDestroy(&CAdiaginv));
82*9fd865deSAlexander       } else {
83*9fd865deSAlexander         if (!lsc->L) {
84*9fd865deSAlexander           PetscCall(MatProductCreate(C, B, NULL, &lsc->L));
85*9fd865deSAlexander           PetscCall(MatProductSetType(lsc->L, MATPRODUCT_AB));
86*9fd865deSAlexander           PetscCall(MatProductSetFromOptions(lsc->L));
87*9fd865deSAlexander           PetscCall(MatProductSymbolic(lsc->L));
88*9fd865deSAlexander         }
89*9fd865deSAlexander         PetscCall(MatProductNumeric(lsc->L));
907e8cb189SBarry Smith       }
917e8cb189SBarry Smith       Lp = L = lsc->L;
927e8cb189SBarry Smith     }
93daa6d063SJed Brown   }
94*9fd865deSAlexander 
959566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(lsc->kspL, L, Lp));
969566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(lsc->kspL));
97*9fd865deSAlexander   if (lsc->commute) {
98*9fd865deSAlexander     PetscCall(KSPSetOperators(lsc->kspMass, Qscale, Qscale));
99*9fd865deSAlexander     PetscCall(KSPSetFromOptions(lsc->kspMass));
100*9fd865deSAlexander   }
1013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
102daa6d063SJed Brown }
103daa6d063SJed Brown 
104d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LSC(PC pc, Vec x, Vec y)
105d71ae5a4SJacob Faibussowitsch {
106daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
107daa6d063SJed Brown   Mat     A, B, C;
108daa6d063SJed Brown 
109daa6d063SJed Brown   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, &B, &C, NULL));
111*9fd865deSAlexander   if (lsc->commute) {
112*9fd865deSAlexander     PetscCall(KSPSolve(lsc->kspMass, x, lsc->Svec0));
113*9fd865deSAlexander     PetscCall(KSPCheckSolve(lsc->kspMass, pc, lsc->Svec0));
114*9fd865deSAlexander     PetscCall(MatMult(B, lsc->Svec0, lsc->Avec0));
115*9fd865deSAlexander     PetscCall(KSPSolve(lsc->kspL, lsc->Avec0, lsc->Avec1));
116*9fd865deSAlexander     PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->Avec1));
117*9fd865deSAlexander     PetscCall(MatMult(A, lsc->Avec1, lsc->Avec0));
118*9fd865deSAlexander     PetscCall(KSPSolve(lsc->kspL, lsc->Avec0, lsc->Avec1));
119*9fd865deSAlexander     PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->Avec1));
120*9fd865deSAlexander     PetscCall(MatMult(C, lsc->Avec1, lsc->Svec0));
121*9fd865deSAlexander     PetscCall(KSPSolve(lsc->kspMass, lsc->Svec0, y));
122*9fd865deSAlexander     PetscCall(KSPCheckSolve(lsc->kspMass, pc, y));
123*9fd865deSAlexander   } else {
124*9fd865deSAlexander     PetscCall(KSPSolve(lsc->kspL, x, lsc->Svec0));
125*9fd865deSAlexander     PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->Svec0));
126*9fd865deSAlexander     PetscCall(MatMult(B, lsc->Svec0, lsc->Avec0));
127*9fd865deSAlexander     if (lsc->scale) PetscCall(VecPointwiseMult(lsc->Avec0, lsc->Avec0, lsc->scale));
128*9fd865deSAlexander     PetscCall(MatMult(A, lsc->Avec0, lsc->Avec1));
129*9fd865deSAlexander     if (lsc->scale) PetscCall(VecPointwiseMult(lsc->Avec1, lsc->Avec1, lsc->scale));
130*9fd865deSAlexander     PetscCall(MatMult(C, lsc->Avec1, lsc->Svec0));
131*9fd865deSAlexander     PetscCall(KSPSolve(lsc->kspL, lsc->Svec0, y));
1329566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(lsc->kspL, pc, y));
133*9fd865deSAlexander   }
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135daa6d063SJed Brown }
136daa6d063SJed Brown 
137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LSC(PC pc)
138d71ae5a4SJacob Faibussowitsch {
139daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
140daa6d063SJed Brown 
141daa6d063SJed Brown   PetscFunctionBegin;
142*9fd865deSAlexander   PetscCall(VecDestroy(&lsc->Avec0));
143*9fd865deSAlexander   PetscCall(VecDestroy(&lsc->Avec1));
144*9fd865deSAlexander   PetscCall(VecDestroy(&lsc->Svec0));
1459566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&lsc->kspL));
146*9fd865deSAlexander   if (lsc->commute) PetscCall(KSPDestroy(&lsc->kspMass));
147*9fd865deSAlexander   if (lsc->L) PetscCall(MatDestroy(&lsc->L));
148*9fd865deSAlexander   if (lsc->scale) PetscCall(VecDestroy(&lsc->scale));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150a06653b4SBarry Smith }
151a06653b4SBarry Smith 
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LSC(PC pc)
153d71ae5a4SJacob Faibussowitsch {
154a06653b4SBarry Smith   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(PCReset_LSC(pc));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
158daa6d063SJed Brown }
159daa6d063SJed Brown 
160d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_LSC(PC pc, PetscOptionItems *PetscOptionsObject)
161d71ae5a4SJacob Faibussowitsch {
162daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
163daa6d063SJed Brown 
164daa6d063SJed Brown   PetscFunctionBegin;
165d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "LSC options");
166d71ae5a4SJacob Faibussowitsch   {
167*9fd865deSAlexander     PetscCall(PetscOptionsBool("-pc_lsc_commute", "Whether to commute the LSC preconditioner in the style of Olshanskii", "None", lsc->commute, &lsc->commute, NULL));
168*9fd865deSAlexander     PetscCall(PetscOptionsBool("-pc_lsc_scale_diag", "Whether to scale BBt products. Will use the inverse of the diagonal of Qscale or A if the former is not provided.", "None", lsc->scalediag, &lsc->scalediag, NULL));
169*9fd865deSAlexander     PetscCheck(!lsc->scalediag || !lsc->commute, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Diagonal-based scaling is not used when doing a commuted LSC. Either do not ask for diagonal-based scaling or use non-commuted LSC in the original style of Elman");
170d71ae5a4SJacob Faibussowitsch   }
171d0609cedSBarry Smith   PetscOptionsHeadEnd();
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173daa6d063SJed Brown }
174daa6d063SJed Brown 
175d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_LSC(PC pc, PetscViewer viewer)
176d71ae5a4SJacob Faibussowitsch {
1777e8cb189SBarry Smith   PC_LSC   *jac = (PC_LSC *)pc->data;
1787e8cb189SBarry Smith   PetscBool iascii;
1797e8cb189SBarry Smith 
1807e8cb189SBarry Smith   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1827e8cb189SBarry Smith   if (iascii) {
1839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
184ca0c957dSBarry Smith     if (jac->kspL) {
1859566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->kspL, viewer));
186ca0c957dSBarry Smith     } else {
1879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC KSP object not yet created, hence cannot display"));
188ca0c957dSBarry Smith     }
189*9fd865deSAlexander     if (jac->commute) {
190*9fd865deSAlexander       if (jac->kspMass) {
191*9fd865deSAlexander         PetscCall(KSPView(jac->kspMass, viewer));
192*9fd865deSAlexander       } else {
193*9fd865deSAlexander         PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC Mass KSP object not yet created, hence cannot display"));
194*9fd865deSAlexander       }
195*9fd865deSAlexander     }
1969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
19711aeaf0aSBarry Smith   }
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1997e8cb189SBarry Smith }
2007e8cb189SBarry Smith 
201daa6d063SJed Brown /*MC
202daa6d063SJed Brown      PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators
203daa6d063SJed Brown 
204daa6d063SJed Brown    Options Database Key:
205daa6d063SJed Brown .    -pc_lsc_scale_diag - Use the diagonal of A for scaling
206daa6d063SJed Brown 
207daa6d063SJed Brown    Level: intermediate
208daa6d063SJed Brown 
209daa6d063SJed Brown    Notes:
210f1580f4eSBarry Smith    This preconditioner will normally be used with `PCFIELDSPLIT` to precondition the Schur complement, but
211daa6d063SJed Brown    it can be used for any Schur complement system.  Consider the Schur complement
212daa6d063SJed Brown 
213daa6d063SJed Brown .vb
2147e8cb189SBarry Smith    S = A11 - A10 inv(A00) A01
215daa6d063SJed Brown .ve
216daa6d063SJed Brown 
217f1580f4eSBarry Smith    `PCLSC` currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
218daa6d063SJed Brown    inv(S) is given by
219daa6d063SJed Brown 
220daa6d063SJed Brown .vb
2217e8cb189SBarry Smith    inv(A10 A01) A10 A00 A01 inv(A10 A01)
222daa6d063SJed Brown .ve
223daa6d063SJed Brown 
2247e8cb189SBarry Smith    The product A10 A01 can be computed for you, but you can provide it (this is
2251c6f1693SPatrick Sanan    usually more efficient anyway).  In the case of incompressible flow, A10 A01 is a Laplacian; call it L.  The current
226daa6d063SJed Brown    interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.
227daa6d063SJed Brown 
228f1580f4eSBarry Smith    If you had called `KSPSetOperators`(ksp,S,Sp), S should have type `MATSCHURCOMPLEMENT` and Sp can be any type you
229f1580f4eSBarry Smith    like (`PCLSC` doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
230daa6d063SJed Brown    For example, you might have setup code like this
231daa6d063SJed Brown 
232daa6d063SJed Brown .vb
233daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
234daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
235daa6d063SJed Brown .ve
236daa6d063SJed Brown 
237daa6d063SJed Brown    And then your Jacobian assembly would look like
238daa6d063SJed Brown 
239daa6d063SJed Brown .vb
240daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
241daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
242daa6d063SJed Brown    if (L) { assembly L }
243daa6d063SJed Brown    if (Lp) { assemble Lp }
244daa6d063SJed Brown .ve
245daa6d063SJed Brown 
246daa6d063SJed Brown    With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L
247daa6d063SJed Brown 
248daa6d063SJed Brown .vb
249daa6d063SJed Brown    -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
250daa6d063SJed Brown .ve
251daa6d063SJed Brown 
252daa6d063SJed Brown    Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.
253daa6d063SJed Brown 
254a674a499SJed Brown    References:
255606c0280SSatish Balay +  * - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006.
256606c0280SSatish Balay -  * - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001.
257a674a499SJed Brown 
258db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `Block_Preconditioners`, `PCFIELDSPLIT`,
259db781477SPatrick Sanan           `PCFieldSplitGetSubKSP()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`,
260f1580f4eSBarry Smith           `MatCreateSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetSubMatrices()`, `MatSchurComplementUpdateSubMatrices()`,
261f1580f4eSBarry Smith           `MatSchurComplementSetAinvType()`, `MatGetSchurComplement()`
262daa6d063SJed Brown M*/
263daa6d063SJed Brown 
264d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc)
265d71ae5a4SJacob Faibussowitsch {
266daa6d063SJed Brown   PC_LSC *lsc;
267daa6d063SJed Brown 
268daa6d063SJed Brown   PetscFunctionBegin;
2694dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lsc));
270daa6d063SJed Brown   pc->data = (void *)lsc;
271daa6d063SJed Brown 
272daa6d063SJed Brown   pc->ops->apply           = PCApply_LSC;
2730a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
274daa6d063SJed Brown   pc->ops->setup           = PCSetUp_LSC;
275a06653b4SBarry Smith   pc->ops->reset           = PCReset_LSC;
276daa6d063SJed Brown   pc->ops->destroy         = PCDestroy_LSC;
277daa6d063SJed Brown   pc->ops->setfromoptions  = PCSetFromOptions_LSC;
2787e8cb189SBarry Smith   pc->ops->view            = PCView_LSC;
2790a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
281daa6d063SJed Brown }
282