xref: /petsc/src/ksp/pc/impls/lsc/lsc.c (revision 3821be0afdf874bf2c399643c829ed821dc5dd56)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
2daa6d063SJed Brown 
3daa6d063SJed Brown typedef struct {
4ace3abfcSBarry Smith   PetscBool allocated;
5ace3abfcSBarry Smith   PetscBool scalediag;
6daa6d063SJed Brown   KSP       kspL;
7daa6d063SJed Brown   Vec       scale;
8daa6d063SJed Brown   Vec       x0, y0, x1;
97e8cb189SBarry Smith   Mat       L; /* keep a copy to reuse when obtained with L = A10*A01 */
10daa6d063SJed Brown } PC_LSC;
11daa6d063SJed Brown 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCLSCAllocate_Private(PC pc)
13d71ae5a4SJacob Faibussowitsch {
14daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
15daa6d063SJed Brown   Mat     A;
16daa6d063SJed Brown 
17daa6d063SJed Brown   PetscFunctionBegin;
183ba16761SJacob Faibussowitsch   if (lsc->allocated) PetscFunctionReturn(PETSC_SUCCESS);
199566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspL));
20*3821be0aSBarry Smith   PetscCall(KSPSetNestLevel(lsc->kspL, pc->kspnestlevel));
219566063dSJacob Faibussowitsch   PetscCall(KSPSetErrorIfNotConverged(lsc->kspL, pc->erroriffailure));
229566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lsc->kspL, (PetscObject)pc, 1));
239566063dSJacob Faibussowitsch   PetscCall(KSPSetType(lsc->kspL, KSPPREONLY));
249566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(lsc->kspL, ((PetscObject)pc)->prefix));
259566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(lsc->kspL, "lsc_"));
269566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, NULL, NULL, NULL));
279566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lsc->x0, &lsc->y0));
289566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &lsc->x1, NULL));
2948a46eb9SPierre Jolivet   if (lsc->scalediag) PetscCall(VecDuplicate(lsc->x0, &lsc->scale));
30daa6d063SJed Brown   lsc->allocated = PETSC_TRUE;
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32daa6d063SJed Brown }
33daa6d063SJed Brown 
34d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LSC(PC pc)
35d71ae5a4SJacob Faibussowitsch {
36daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
377e8cb189SBarry Smith   Mat     L, Lp, B, C;
38daa6d063SJed Brown 
39daa6d063SJed Brown   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PCLSCAllocate_Private(pc));
419566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_L", (PetscObject *)&L));
429566063dSJacob Faibussowitsch   if (!L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_L", (PetscObject *)&L));
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Lp", (PetscObject *)&Lp));
449566063dSJacob Faibussowitsch   if (!Lp) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Lp", (PetscObject *)&Lp));
457e8cb189SBarry Smith   if (!L) {
469566063dSJacob Faibussowitsch     PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, NULL, &B, &C, NULL));
477e8cb189SBarry Smith     if (!lsc->L) {
489566063dSJacob Faibussowitsch       PetscCall(MatMatMult(C, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &lsc->L));
497e8cb189SBarry Smith     } else {
509566063dSJacob Faibussowitsch       PetscCall(MatMatMult(C, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &lsc->L));
517e8cb189SBarry Smith     }
527e8cb189SBarry Smith     Lp = L = lsc->L;
537e8cb189SBarry Smith   }
54daa6d063SJed Brown   if (lsc->scale) {
55daa6d063SJed Brown     Mat Ap;
569566063dSJacob Faibussowitsch     PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, &Ap, NULL, NULL, NULL));
579566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Ap, lsc->scale)); /* Should be the mass matrix, but we don't have plumbing for that yet */
589566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(lsc->scale));
59daa6d063SJed Brown   }
609566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(lsc->kspL, L, Lp));
619566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(lsc->kspL));
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63daa6d063SJed Brown }
64daa6d063SJed Brown 
65d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LSC(PC pc, Vec x, Vec y)
66d71ae5a4SJacob Faibussowitsch {
67daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
68daa6d063SJed Brown   Mat     A, B, C;
69daa6d063SJed Brown 
70daa6d063SJed Brown   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, &B, &C, NULL));
729566063dSJacob Faibussowitsch   PetscCall(KSPSolve(lsc->kspL, x, lsc->x1));
739566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->x1));
749566063dSJacob Faibussowitsch   PetscCall(MatMult(B, lsc->x1, lsc->x0));
751baa6e33SBarry Smith   if (lsc->scale) PetscCall(VecPointwiseMult(lsc->x0, lsc->x0, lsc->scale));
769566063dSJacob Faibussowitsch   PetscCall(MatMult(A, lsc->x0, lsc->y0));
771baa6e33SBarry Smith   if (lsc->scale) PetscCall(VecPointwiseMult(lsc->y0, lsc->y0, lsc->scale));
789566063dSJacob Faibussowitsch   PetscCall(MatMult(C, lsc->y0, lsc->x1));
799566063dSJacob Faibussowitsch   PetscCall(KSPSolve(lsc->kspL, lsc->x1, y));
809566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(lsc->kspL, pc, y));
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82daa6d063SJed Brown }
83daa6d063SJed Brown 
84d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LSC(PC pc)
85d71ae5a4SJacob Faibussowitsch {
86daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
87daa6d063SJed Brown 
88daa6d063SJed Brown   PetscFunctionBegin;
899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->x0));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->y0));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->x1));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->scale));
939566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&lsc->kspL));
949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lsc->L));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96a06653b4SBarry Smith }
97a06653b4SBarry Smith 
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LSC(PC pc)
99d71ae5a4SJacob Faibussowitsch {
100a06653b4SBarry Smith   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(PCReset_LSC(pc));
1029566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104daa6d063SJed Brown }
105daa6d063SJed Brown 
106d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_LSC(PC pc, PetscOptionItems *PetscOptionsObject)
107d71ae5a4SJacob Faibussowitsch {
108daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
109daa6d063SJed Brown 
110daa6d063SJed Brown   PetscFunctionBegin;
111d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "LSC options");
112d71ae5a4SJacob Faibussowitsch   {
113d71ae5a4SJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_lsc_scale_diag", "Use diagonal of velocity block (A) for scaling", "None", lsc->scalediag, &lsc->scalediag, NULL));
114d71ae5a4SJacob Faibussowitsch   }
115d0609cedSBarry Smith   PetscOptionsHeadEnd();
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117daa6d063SJed Brown }
118daa6d063SJed Brown 
119d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_LSC(PC pc, PetscViewer viewer)
120d71ae5a4SJacob Faibussowitsch {
1217e8cb189SBarry Smith   PC_LSC   *jac = (PC_LSC *)pc->data;
1227e8cb189SBarry Smith   PetscBool iascii;
1237e8cb189SBarry Smith 
1247e8cb189SBarry Smith   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1267e8cb189SBarry Smith   if (iascii) {
1279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
128ca0c957dSBarry Smith     if (jac->kspL) {
1299566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->kspL, viewer));
130ca0c957dSBarry Smith     } else {
1319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC KSP object not yet created, hence cannot display"));
132ca0c957dSBarry Smith     }
1339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
13411aeaf0aSBarry Smith   }
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1367e8cb189SBarry Smith }
1377e8cb189SBarry Smith 
138daa6d063SJed Brown /*MC
139daa6d063SJed Brown      PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators
140daa6d063SJed Brown 
141daa6d063SJed Brown    Options Database Key:
142daa6d063SJed Brown .    -pc_lsc_scale_diag - Use the diagonal of A for scaling
143daa6d063SJed Brown 
144daa6d063SJed Brown    Level: intermediate
145daa6d063SJed Brown 
146daa6d063SJed Brown    Notes:
147f1580f4eSBarry Smith    This preconditioner will normally be used with `PCFIELDSPLIT` to precondition the Schur complement, but
148daa6d063SJed Brown    it can be used for any Schur complement system.  Consider the Schur complement
149daa6d063SJed Brown 
150daa6d063SJed Brown .vb
1517e8cb189SBarry Smith    S = A11 - A10 inv(A00) A01
152daa6d063SJed Brown .ve
153daa6d063SJed Brown 
154f1580f4eSBarry Smith    `PCLSC` currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
155daa6d063SJed Brown    inv(S) is given by
156daa6d063SJed Brown 
157daa6d063SJed Brown .vb
1587e8cb189SBarry Smith    inv(A10 A01) A10 A00 A01 inv(A10 A01)
159daa6d063SJed Brown .ve
160daa6d063SJed Brown 
1617e8cb189SBarry Smith    The product A10 A01 can be computed for you, but you can provide it (this is
1621c6f1693SPatrick Sanan    usually more efficient anyway).  In the case of incompressible flow, A10 A01 is a Laplacian; call it L.  The current
163daa6d063SJed Brown    interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.
164daa6d063SJed Brown 
165f1580f4eSBarry Smith    If you had called `KSPSetOperators`(ksp,S,Sp), S should have type `MATSCHURCOMPLEMENT` and Sp can be any type you
166f1580f4eSBarry Smith    like (`PCLSC` doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
167daa6d063SJed Brown    For example, you might have setup code like this
168daa6d063SJed Brown 
169daa6d063SJed Brown .vb
170daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
171daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
172daa6d063SJed Brown .ve
173daa6d063SJed Brown 
174daa6d063SJed Brown    And then your Jacobian assembly would look like
175daa6d063SJed Brown 
176daa6d063SJed Brown .vb
177daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
178daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
179daa6d063SJed Brown    if (L) { assembly L }
180daa6d063SJed Brown    if (Lp) { assemble Lp }
181daa6d063SJed Brown .ve
182daa6d063SJed Brown 
183daa6d063SJed Brown    With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L
184daa6d063SJed Brown 
185daa6d063SJed Brown .vb
186daa6d063SJed Brown    -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
187daa6d063SJed Brown .ve
188daa6d063SJed Brown 
189daa6d063SJed Brown    Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.
190daa6d063SJed Brown 
191a674a499SJed Brown    References:
192606c0280SSatish Balay +  * - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006.
193606c0280SSatish Balay -  * - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001.
194a674a499SJed Brown 
195db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `Block_Preconditioners`, `PCFIELDSPLIT`,
196db781477SPatrick Sanan           `PCFieldSplitGetSubKSP()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`,
197f1580f4eSBarry Smith           `MatCreateSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetSubMatrices()`, `MatSchurComplementUpdateSubMatrices()`,
198f1580f4eSBarry Smith           `MatSchurComplementSetAinvType()`, `MatGetSchurComplement()`
199daa6d063SJed Brown M*/
200daa6d063SJed Brown 
201d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc)
202d71ae5a4SJacob Faibussowitsch {
203daa6d063SJed Brown   PC_LSC *lsc;
204daa6d063SJed Brown 
205daa6d063SJed Brown   PetscFunctionBegin;
2064dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lsc));
207daa6d063SJed Brown   pc->data = (void *)lsc;
208daa6d063SJed Brown 
209daa6d063SJed Brown   pc->ops->apply           = PCApply_LSC;
2100a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
211daa6d063SJed Brown   pc->ops->setup           = PCSetUp_LSC;
212a06653b4SBarry Smith   pc->ops->reset           = PCReset_LSC;
213daa6d063SJed Brown   pc->ops->destroy         = PCDestroy_LSC;
214daa6d063SJed Brown   pc->ops->setfromoptions  = PCSetFromOptions_LSC;
2157e8cb189SBarry Smith   pc->ops->view            = PCView_LSC;
2160a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218daa6d063SJed Brown }
219