xref: /petsc/src/ksp/pc/impls/lsc/lsc.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
12daa6d063SJed Brown static PetscErrorCode PCLSCAllocate_Private(PC pc)
13daa6d063SJed Brown {
14daa6d063SJed Brown   PC_LSC         *lsc = (PC_LSC*)pc->data;
15daa6d063SJed Brown   Mat            A;
16daa6d063SJed Brown 
17daa6d063SJed Brown   PetscFunctionBegin;
18daa6d063SJed Brown   if (lsc->allocated) PetscFunctionReturn(0);
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));
28daa6d063SJed Brown   if (lsc->scalediag) {
299566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(lsc->x0,&lsc->scale));
30daa6d063SJed Brown   }
31daa6d063SJed Brown   lsc->allocated = PETSC_TRUE;
32daa6d063SJed Brown   PetscFunctionReturn(0);
33daa6d063SJed Brown }
34daa6d063SJed Brown 
35daa6d063SJed Brown static PetscErrorCode PCSetUp_LSC(PC pc)
36daa6d063SJed Brown {
37daa6d063SJed Brown   PC_LSC         *lsc = (PC_LSC*)pc->data;
387e8cb189SBarry Smith   Mat            L,Lp,B,C;
39daa6d063SJed Brown 
40daa6d063SJed Brown   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(PCLSCAllocate_Private(pc));
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->mat,"LSC_L",(PetscObject*)&L));
439566063dSJacob Faibussowitsch   if (!L) PetscCall(PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L));
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp));
459566063dSJacob Faibussowitsch   if (!Lp) PetscCall(PetscObjectQuery((PetscObject)pc->mat,"LSC_Lp",(PetscObject*)&Lp));
467e8cb189SBarry Smith   if (!L) {
479566063dSJacob Faibussowitsch     PetscCall(MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL));
487e8cb189SBarry Smith     if (!lsc->L) {
499566063dSJacob Faibussowitsch       PetscCall(MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L));
507e8cb189SBarry Smith     } else {
519566063dSJacob Faibussowitsch       PetscCall(MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L));
527e8cb189SBarry Smith     }
537e8cb189SBarry Smith     Lp = L = lsc->L;
547e8cb189SBarry Smith   }
55daa6d063SJed Brown   if (lsc->scale) {
56daa6d063SJed Brown     Mat Ap;
579566063dSJacob Faibussowitsch     PetscCall(MatSchurComplementGetSubMatrices(pc->mat,NULL,&Ap,NULL,NULL,NULL));
589566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Ap,lsc->scale)); /* Should be the mass matrix, but we don't have plumbing for that yet */
599566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(lsc->scale));
60daa6d063SJed Brown   }
619566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(lsc->kspL,L,Lp));
629566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(lsc->kspL));
63daa6d063SJed Brown   PetscFunctionReturn(0);
64daa6d063SJed Brown }
65daa6d063SJed Brown 
66daa6d063SJed Brown static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y)
67daa6d063SJed Brown {
68daa6d063SJed Brown   PC_LSC         *lsc = (PC_LSC*)pc->data;
69daa6d063SJed Brown   Mat            A,B,C;
70daa6d063SJed Brown 
71daa6d063SJed Brown   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(MatSchurComplementGetSubMatrices(pc->mat,&A,NULL,&B,&C,NULL));
739566063dSJacob Faibussowitsch   PetscCall(KSPSolve(lsc->kspL,x,lsc->x1));
749566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(lsc->kspL,pc,lsc->x1));
759566063dSJacob Faibussowitsch   PetscCall(MatMult(B,lsc->x1,lsc->x0));
76daa6d063SJed Brown   if (lsc->scale) {
779566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale));
78daa6d063SJed Brown   }
799566063dSJacob Faibussowitsch   PetscCall(MatMult(A,lsc->x0,lsc->y0));
80daa6d063SJed Brown   if (lsc->scale) {
819566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale));
82daa6d063SJed Brown   }
839566063dSJacob Faibussowitsch   PetscCall(MatMult(C,lsc->y0,lsc->x1));
849566063dSJacob Faibussowitsch   PetscCall(KSPSolve(lsc->kspL,lsc->x1,y));
859566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(lsc->kspL,pc,y));
86daa6d063SJed Brown   PetscFunctionReturn(0);
87daa6d063SJed Brown }
88daa6d063SJed Brown 
89a06653b4SBarry Smith static PetscErrorCode PCReset_LSC(PC pc)
90daa6d063SJed Brown {
91daa6d063SJed Brown   PC_LSC         *lsc = (PC_LSC*)pc->data;
92daa6d063SJed Brown 
93daa6d063SJed Brown   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->x0));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->y0));
969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->x1));
979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lsc->scale));
989566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&lsc->kspL));
999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lsc->L));
100a06653b4SBarry Smith   PetscFunctionReturn(0);
101a06653b4SBarry Smith }
102a06653b4SBarry Smith 
103a06653b4SBarry Smith static PetscErrorCode PCDestroy_LSC(PC pc)
104a06653b4SBarry Smith {
105a06653b4SBarry Smith   PetscFunctionBegin;
1069566063dSJacob Faibussowitsch   PetscCall(PCReset_LSC(pc));
1079566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
108daa6d063SJed Brown   PetscFunctionReturn(0);
109daa6d063SJed Brown }
110daa6d063SJed Brown 
1114416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LSC(PetscOptionItems *PetscOptionsObject,PC pc)
112daa6d063SJed Brown {
113daa6d063SJed Brown   PC_LSC         *lsc = (PC_LSC*)pc->data;
114daa6d063SJed Brown 
115daa6d063SJed Brown   PetscFunctionBegin;
116*d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"LSC options");
117e55864a3SBarry Smith   {
1189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,NULL));
1191a1499c8SBarry Smith   }
120*d0609cedSBarry Smith   PetscOptionsHeadEnd();
121daa6d063SJed Brown   PetscFunctionReturn(0);
122daa6d063SJed Brown }
123daa6d063SJed Brown 
1247e8cb189SBarry Smith static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer)
1257e8cb189SBarry Smith {
1267e8cb189SBarry Smith   PC_LSC         *jac = (PC_LSC*)pc->data;
1277e8cb189SBarry Smith   PetscBool      iascii;
1287e8cb189SBarry Smith 
1297e8cb189SBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1317e8cb189SBarry Smith   if (iascii) {
1329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
133ca0c957dSBarry Smith     if (jac->kspL) {
1349566063dSJacob Faibussowitsch       PetscCall(KSPView(jac->kspL,viewer));
135ca0c957dSBarry Smith     } else {
1369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"PCLSC KSP object not yet created, hence cannot display"));
137ca0c957dSBarry Smith     }
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
13911aeaf0aSBarry Smith   }
1407e8cb189SBarry Smith   PetscFunctionReturn(0);
1417e8cb189SBarry Smith }
1427e8cb189SBarry Smith 
143daa6d063SJed Brown /*MC
144daa6d063SJed Brown      PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators
145daa6d063SJed Brown 
146daa6d063SJed Brown    Options Database Key:
147daa6d063SJed Brown .    -pc_lsc_scale_diag - Use the diagonal of A for scaling
148daa6d063SJed Brown 
149daa6d063SJed Brown    Level: intermediate
150daa6d063SJed Brown 
151daa6d063SJed Brown    Notes:
152daa6d063SJed Brown    This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but
153daa6d063SJed Brown    it can be used for any Schur complement system.  Consider the Schur complement
154daa6d063SJed Brown 
155daa6d063SJed Brown .vb
1567e8cb189SBarry Smith    S = A11 - A10 inv(A00) A01
157daa6d063SJed Brown .ve
158daa6d063SJed Brown 
1597e8cb189SBarry Smith    PCLSC currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
160daa6d063SJed Brown    inv(S) is given by
161daa6d063SJed Brown 
162daa6d063SJed Brown .vb
1637e8cb189SBarry Smith    inv(A10 A01) A10 A00 A01 inv(A10 A01)
164daa6d063SJed Brown .ve
165daa6d063SJed Brown 
1667e8cb189SBarry Smith    The product A10 A01 can be computed for you, but you can provide it (this is
1671c6f1693SPatrick Sanan    usually more efficient anyway).  In the case of incompressible flow, A10 A01 is a Laplacian; call it L.  The current
168daa6d063SJed Brown    interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.
169daa6d063SJed Brown 
17062279869SBarry Smith    If you had called KSPSetOperators(ksp,S,Sp), S should have type MATSCHURCOMPLEMENT and Sp can be any type you
171daa6d063SJed Brown    like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
172daa6d063SJed Brown    For example, you might have setup code like this
173daa6d063SJed Brown 
174daa6d063SJed Brown .vb
175daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
176daa6d063SJed Brown    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
177daa6d063SJed Brown .ve
178daa6d063SJed Brown 
179daa6d063SJed Brown    And then your Jacobian assembly would look like
180daa6d063SJed Brown 
181daa6d063SJed Brown .vb
182daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
183daa6d063SJed Brown    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
184daa6d063SJed Brown    if (L) { assembly L }
185daa6d063SJed Brown    if (Lp) { assemble Lp }
186daa6d063SJed Brown .ve
187daa6d063SJed Brown 
188daa6d063SJed Brown    With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L
189daa6d063SJed Brown 
190daa6d063SJed Brown .vb
191daa6d063SJed Brown    -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
192daa6d063SJed Brown .ve
193daa6d063SJed Brown 
194daa6d063SJed Brown    Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.
195daa6d063SJed Brown 
196a674a499SJed Brown    References:
197606c0280SSatish Balay +  * - Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006.
198606c0280SSatish Balay -  * - Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier Stokes equations for incompressible flow, 2001.
199a674a499SJed Brown 
200daa6d063SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT,
20129f8a81cSJed Brown            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
202daa6d063SJed Brown            MatCreateSchurComplement()
203daa6d063SJed Brown M*/
204daa6d063SJed Brown 
2058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LSC(PC pc)
206daa6d063SJed Brown {
207daa6d063SJed Brown   PC_LSC         *lsc;
208daa6d063SJed Brown 
209daa6d063SJed Brown   PetscFunctionBegin;
2109566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&lsc));
211daa6d063SJed Brown   pc->data = (void*)lsc;
212daa6d063SJed Brown 
213daa6d063SJed Brown   pc->ops->apply           = PCApply_LSC;
2140a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
215daa6d063SJed Brown   pc->ops->setup           = PCSetUp_LSC;
216a06653b4SBarry Smith   pc->ops->reset           = PCReset_LSC;
217daa6d063SJed Brown   pc->ops->destroy         = PCDestroy_LSC;
218daa6d063SJed Brown   pc->ops->setfromoptions  = PCSetFromOptions_LSC;
2197e8cb189SBarry Smith   pc->ops->view            = PCView_LSC;
2200a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
221daa6d063SJed Brown   PetscFunctionReturn(0);
222daa6d063SJed Brown }
223