xref: /petsc/src/ksp/pc/impls/lsc/lsc.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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;
18*3ba16761SJacob Faibussowitsch   if (lsc->allocated) PetscFunctionReturn(PETSC_SUCCESS);
199566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &lsc->kspL));
209566063dSJacob Faibussowitsch   PetscCall(KSPSetErrorIfNotConverged(lsc->kspL, pc->erroriffailure));
219566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lsc->kspL, (PetscObject)pc, 1));
229566063dSJacob Faibussowitsch   PetscCall(KSPSetType(lsc->kspL, KSPPREONLY));
239566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(lsc->kspL, ((PetscObject)pc)->prefix));
249566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(lsc->kspL, "lsc_"));
259566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, NULL, NULL, NULL));
269566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lsc->x0, &lsc->y0));
279566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &lsc->x1, NULL));
2848a46eb9SPierre Jolivet   if (lsc->scalediag) PetscCall(VecDuplicate(lsc->x0, &lsc->scale));
29daa6d063SJed Brown   lsc->allocated = PETSC_TRUE;
30*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31daa6d063SJed Brown }
32daa6d063SJed Brown 
33d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LSC(PC pc)
34d71ae5a4SJacob Faibussowitsch {
35daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
367e8cb189SBarry Smith   Mat     L, Lp, B, C;
37daa6d063SJed Brown 
38daa6d063SJed Brown   PetscFunctionBegin;
399566063dSJacob Faibussowitsch   PetscCall(PCLSCAllocate_Private(pc));
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_L", (PetscObject *)&L));
419566063dSJacob Faibussowitsch   if (!L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_L", (PetscObject *)&L));
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->pmat, "LSC_Lp", (PetscObject *)&Lp));
439566063dSJacob Faibussowitsch   if (!Lp) PetscCall(PetscObjectQuery((PetscObject)pc->mat, "LSC_Lp", (PetscObject *)&Lp));
447e8cb189SBarry Smith   if (!L) {
459566063dSJacob Faibussowitsch     PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, NULL, &B, &C, NULL));
467e8cb189SBarry Smith     if (!lsc->L) {
479566063dSJacob Faibussowitsch       PetscCall(MatMatMult(C, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &lsc->L));
487e8cb189SBarry Smith     } else {
499566063dSJacob Faibussowitsch       PetscCall(MatMatMult(C, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &lsc->L));
507e8cb189SBarry Smith     }
517e8cb189SBarry Smith     Lp = L = lsc->L;
527e8cb189SBarry Smith   }
53daa6d063SJed Brown   if (lsc->scale) {
54daa6d063SJed Brown     Mat Ap;
559566063dSJacob Faibussowitsch     PetscCall(MatSchurComplementGetSubMatrices(pc->mat, NULL, &Ap, NULL, NULL, NULL));
569566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Ap, lsc->scale)); /* Should be the mass matrix, but we don't have plumbing for that yet */
579566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(lsc->scale));
58daa6d063SJed Brown   }
599566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(lsc->kspL, L, Lp));
609566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(lsc->kspL));
61*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62daa6d063SJed Brown }
63daa6d063SJed Brown 
64d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LSC(PC pc, Vec x, Vec y)
65d71ae5a4SJacob Faibussowitsch {
66daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
67daa6d063SJed Brown   Mat     A, B, C;
68daa6d063SJed Brown 
69daa6d063SJed Brown   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat, &A, NULL, &B, &C, NULL));
719566063dSJacob Faibussowitsch   PetscCall(KSPSolve(lsc->kspL, x, lsc->x1));
729566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(lsc->kspL, pc, lsc->x1));
739566063dSJacob Faibussowitsch   PetscCall(MatMult(B, lsc->x1, lsc->x0));
741baa6e33SBarry Smith   if (lsc->scale) PetscCall(VecPointwiseMult(lsc->x0, lsc->x0, lsc->scale));
759566063dSJacob Faibussowitsch   PetscCall(MatMult(A, lsc->x0, lsc->y0));
761baa6e33SBarry Smith   if (lsc->scale) PetscCall(VecPointwiseMult(lsc->y0, lsc->y0, lsc->scale));
779566063dSJacob Faibussowitsch   PetscCall(MatMult(C, lsc->y0, lsc->x1));
789566063dSJacob Faibussowitsch   PetscCall(KSPSolve(lsc->kspL, lsc->x1, y));
799566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(lsc->kspL, pc, y));
80*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81daa6d063SJed Brown }
82daa6d063SJed Brown 
83d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LSC(PC pc)
84d71ae5a4SJacob Faibussowitsch {
85daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
86daa6d063SJed Brown 
87daa6d063SJed Brown   PetscFunctionBegin;
889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->x0));
899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->y0));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->x1));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->scale));
929566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&lsc->kspL));
939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lsc->L));
94*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95a06653b4SBarry Smith }
96a06653b4SBarry Smith 
97d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LSC(PC pc)
98d71ae5a4SJacob Faibussowitsch {
99a06653b4SBarry Smith   PetscFunctionBegin;
1009566063dSJacob Faibussowitsch   PetscCall(PCReset_LSC(pc));
1019566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
102*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103daa6d063SJed Brown }
104daa6d063SJed Brown 
105d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_LSC(PC pc, PetscOptionItems *PetscOptionsObject)
106d71ae5a4SJacob Faibussowitsch {
107daa6d063SJed Brown   PC_LSC *lsc = (PC_LSC *)pc->data;
108daa6d063SJed Brown 
109daa6d063SJed Brown   PetscFunctionBegin;
110d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "LSC options");
111d71ae5a4SJacob Faibussowitsch   {
112d71ae5a4SJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_lsc_scale_diag", "Use diagonal of velocity block (A) for scaling", "None", lsc->scalediag, &lsc->scalediag, NULL));
113d71ae5a4SJacob Faibussowitsch   }
114d0609cedSBarry Smith   PetscOptionsHeadEnd();
115*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116daa6d063SJed Brown }
117daa6d063SJed Brown 
118d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_LSC(PC pc, PetscViewer viewer)
119d71ae5a4SJacob Faibussowitsch {
1207e8cb189SBarry Smith   PC_LSC   *jac = (PC_LSC *)pc->data;
1217e8cb189SBarry Smith   PetscBool iascii;
1227e8cb189SBarry Smith 
1237e8cb189SBarry Smith   PetscFunctionBegin;
1249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1257e8cb189SBarry Smith   if (iascii) {
1269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
127ca0c957dSBarry Smith     if (jac->kspL) {
1289566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->kspL, viewer));
129ca0c957dSBarry Smith     } else {
1309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "PCLSC KSP object not yet created, hence cannot display"));
131ca0c957dSBarry Smith     }
1329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
13311aeaf0aSBarry Smith   }
134*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1357e8cb189SBarry Smith }
1367e8cb189SBarry Smith 
137daa6d063SJed Brown /*MC
138daa6d063SJed Brown      PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators
139daa6d063SJed Brown 
140daa6d063SJed Brown    Options Database Key:
141daa6d063SJed Brown .    -pc_lsc_scale_diag - Use the diagonal of A for scaling
142daa6d063SJed Brown 
143daa6d063SJed Brown    Level: intermediate
144daa6d063SJed Brown 
145daa6d063SJed Brown    Notes:
146f1580f4eSBarry Smith    This preconditioner will normally be used with `PCFIELDSPLIT` to precondition the Schur complement, but
147daa6d063SJed Brown    it can be used for any Schur complement system.  Consider the Schur complement
148daa6d063SJed Brown 
149daa6d063SJed Brown .vb
1507e8cb189SBarry Smith    S = A11 - A10 inv(A00) A01
151daa6d063SJed Brown .ve
152daa6d063SJed Brown 
153f1580f4eSBarry Smith    `PCLSC` currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
154daa6d063SJed Brown    inv(S) is given by
155daa6d063SJed Brown 
156daa6d063SJed Brown .vb
1577e8cb189SBarry Smith    inv(A10 A01) A10 A00 A01 inv(A10 A01)
158daa6d063SJed Brown .ve
159daa6d063SJed Brown 
1607e8cb189SBarry Smith    The product A10 A01 can be computed for you, but you can provide it (this is
1611c6f1693SPatrick Sanan    usually more efficient anyway).  In the case of incompressible flow, A10 A01 is a Laplacian; call it L.  The current
162daa6d063SJed Brown    interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.
163daa6d063SJed Brown 
164f1580f4eSBarry Smith    If you had called `KSPSetOperators`(ksp,S,Sp), S should have type `MATSCHURCOMPLEMENT` and Sp can be any type you
165f1580f4eSBarry Smith    like (`PCLSC` doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
166daa6d063SJed Brown    For example, you might have setup code like this
167daa6d063SJed Brown 
168daa6d063SJed Brown .vb
169daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
170daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
171daa6d063SJed Brown .ve
172daa6d063SJed Brown 
173daa6d063SJed Brown    And then your Jacobian assembly would look like
174daa6d063SJed Brown 
175daa6d063SJed Brown .vb
176daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
177daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
178daa6d063SJed Brown    if (L) { assembly L }
179daa6d063SJed Brown    if (Lp) { assemble Lp }
180daa6d063SJed Brown .ve
181daa6d063SJed Brown 
182daa6d063SJed Brown    With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L
183daa6d063SJed Brown 
184daa6d063SJed Brown .vb
185daa6d063SJed Brown    -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
186daa6d063SJed Brown .ve
187daa6d063SJed Brown 
188daa6d063SJed Brown    Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.
189daa6d063SJed Brown 
190a674a499SJed Brown    References:
191606c0280SSatish Balay +  * - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006.
192606c0280SSatish Balay -  * - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001.
193a674a499SJed Brown 
194db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `Block_Preconditioners`, `PCFIELDSPLIT`,
195db781477SPatrick Sanan           `PCFieldSplitGetSubKSP()`, `PCFieldSplitSetFields()`, `PCFieldSplitSetType()`, `PCFieldSplitSetIS()`, `PCFieldSplitSetSchurPre()`,
196f1580f4eSBarry Smith           `MatCreateSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetSubMatrices()`, `MatSchurComplementUpdateSubMatrices()`,
197f1580f4eSBarry Smith           `MatSchurComplementSetAinvType()`, `MatGetSchurComplement()`
198daa6d063SJed Brown M*/
199daa6d063SJed Brown 
200d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc)
201d71ae5a4SJacob Faibussowitsch {
202daa6d063SJed Brown   PC_LSC *lsc;
203daa6d063SJed Brown 
204daa6d063SJed Brown   PetscFunctionBegin;
2054dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lsc));
206daa6d063SJed Brown   pc->data = (void *)lsc;
207daa6d063SJed Brown 
208daa6d063SJed Brown   pc->ops->apply           = PCApply_LSC;
2090a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
210daa6d063SJed Brown   pc->ops->setup           = PCSetUp_LSC;
211a06653b4SBarry Smith   pc->ops->reset           = PCReset_LSC;
212daa6d063SJed Brown   pc->ops->destroy         = PCDestroy_LSC;
213daa6d063SJed Brown   pc->ops->setfromoptions  = PCSetFromOptions_LSC;
2147e8cb189SBarry Smith   pc->ops->view            = PCView_LSC;
2150a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
216*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
217daa6d063SJed Brown }
218