xref: /petsc/src/ksp/pc/impls/lsc/lsc.c (revision 01a79839fc82a7dabb7a87cd2a8bb532c6bfa88d)
1 #define PETSCKSP_DLL
2 
3 
4 #include "private/pcimpl.h"   /*I "petscpc.h" I*/
5 
6 typedef struct {
7   PetscBool  allocated;
8   PetscBool  scalediag;
9   KSP        kspL;
10   Vec        scale;
11   Vec        x0,y0,x1;
12   Mat        L;            /* keep a copy to reuse when obtained with L = A10*A01 */
13 } PC_LSC;
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "PCLSCAllocate_Private"
17 static PetscErrorCode PCLSCAllocate_Private(PC pc)
18 {
19   PC_LSC         *lsc = (PC_LSC*)pc->data;
20   Mat             A;
21   PetscErrorCode  ierr;
22 
23   PetscFunctionBegin;
24   if (lsc->allocated) PetscFunctionReturn(0);
25   ierr = KSPCreate(((PetscObject)pc)->comm,&lsc->kspL);CHKERRQ(ierr);
26   ierr = PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);CHKERRQ(ierr);
27   ierr = KSPSetType(lsc->kspL,KSPPREONLY);CHKERRQ(ierr);
28   ierr = KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);CHKERRQ(ierr);
29   ierr = KSPAppendOptionsPrefix(lsc->kspL,"lsc_");CHKERRQ(ierr);
30   ierr = KSPSetFromOptions(lsc->kspL);CHKERRQ(ierr);
31   ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
32   ierr = MatGetVecs(A,&lsc->x0,&lsc->y0);CHKERRQ(ierr);
33   ierr = MatGetVecs(pc->pmat,&lsc->x1,PETSC_NULL);CHKERRQ(ierr);
34   if (lsc->scalediag) {
35     ierr = VecDuplicate(lsc->x0,&lsc->scale);CHKERRQ(ierr);
36   }
37   lsc->allocated = PETSC_TRUE;
38   PetscFunctionReturn(0);
39 }
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "PCSetUp_LSC"
43 static PetscErrorCode PCSetUp_LSC(PC pc)
44 {
45   PC_LSC         *lsc = (PC_LSC*)pc->data;
46   Mat             L,Lp,B,C;
47   PetscErrorCode  ierr;
48 
49   PetscFunctionBegin;
50   ierr = PCLSCAllocate_Private(pc);CHKERRQ(ierr);
51   ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);CHKERRQ(ierr);
52   ierr = PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);CHKERRQ(ierr);
53   if (!L) {
54     ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr);
55     if (!lsc->L) {
56       ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr);
57     } else {
58       ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr);
59     }
60     Lp = L = lsc->L;
61   }
62   if (lsc->scale) {
63     Mat Ap;
64     ierr = MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,&Ap,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
65     ierr = MatGetDiagonal(Ap,lsc->scale);CHKERRQ(ierr); /* Should be the mass matrix, but we don't have plumbing for that yet */
66     ierr = VecReciprocal(lsc->scale);CHKERRQ(ierr);
67   }
68   ierr = KSPSetOperators(lsc->kspL,L,Lp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
69   PetscFunctionReturn(0);
70 }
71 
72 #undef __FUNCT__
73 #define __FUNCT__ "PCApply_LSC"
74 static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y)
75 {
76   PC_LSC        *lsc = (PC_LSC*)pc->data;
77   Mat            A,B,C;
78   PetscErrorCode ierr;
79 
80   PetscFunctionBegin;
81   ierr = MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,&B,&C,PETSC_NULL);CHKERRQ(ierr);
82   ierr = KSPSolve(lsc->kspL,x,lsc->x1);CHKERRQ(ierr);
83   ierr = MatMult(B,lsc->x1,lsc->x0);CHKERRQ(ierr);
84   if (lsc->scale) {
85     ierr = VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);CHKERRQ(ierr);
86   }
87   ierr = MatMult(A,lsc->x0,lsc->y0);CHKERRQ(ierr);
88   if (lsc->scale) {
89     ierr = VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);CHKERRQ(ierr);
90   }
91   ierr = MatMult(C,lsc->y0,lsc->x1);CHKERRQ(ierr);
92   ierr = KSPSolve(lsc->kspL,lsc->x1,y);CHKERRQ(ierr);
93   PetscFunctionReturn(0);
94 }
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "PCDestroy_LSC"
98 static PetscErrorCode PCDestroy_LSC(PC pc)
99 {
100   PC_LSC         *lsc = (PC_LSC*)pc->data;
101   PetscErrorCode ierr;
102 
103   PetscFunctionBegin;
104   if (lsc->x0)    {ierr = VecDestroy(lsc->x0);CHKERRQ(ierr);}
105   if (lsc->y0)    {ierr = VecDestroy(lsc->y0);CHKERRQ(ierr);}
106   if (lsc->x1)    {ierr = VecDestroy(lsc->x1);CHKERRQ(ierr);}
107   if (lsc->scale) {ierr = VecDestroy(lsc->scale);CHKERRQ(ierr);}
108   if (lsc->kspL)  {ierr = KSPDestroy(lsc->kspL);CHKERRQ(ierr);}
109   if (lsc->L)     {ierr = MatDestroy(lsc->L);CHKERRQ(ierr);}
110   ierr = PetscFree(lsc);CHKERRQ(ierr);
111   PetscFunctionReturn(0);
112 }
113 
114 #undef __FUNCT__
115 #define __FUNCT__ "PCSetFromOptions_LSC"
116 static PetscErrorCode PCSetFromOptions_LSC(PC pc)
117 {
118   PC_LSC         *lsc = (PC_LSC*)pc->data;
119   PetscErrorCode  ierr;
120 
121   PetscFunctionBegin;
122   ierr = PetscOptionsHead("LSC options");CHKERRQ(ierr);
123   {
124     ierr = PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,PETSC_NULL);CHKERRQ(ierr);
125   }
126   ierr = PetscOptionsTail();CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "PCView_LSC"
132 static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer)
133 {
134   PC_LSC           *jac = (PC_LSC*)pc->data;
135   PetscErrorCode   ierr;
136   PetscBool        iascii;
137 
138   PetscFunctionBegin;
139   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
140   if (iascii) {
141     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
142     ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr);
143     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
144   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for LSC",((PetscObject)viewer)->type_name);
145   PetscFunctionReturn(0);
146 }
147 
148 /*MC
149      PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators
150 
151    Options Database Key:
152 .    -pc_lsc_scale_diag - Use the diagonal of A for scaling
153 
154    Level: intermediate
155 
156    Notes:
157    This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but
158    it can be used for any Schur complement system.  Consider the Schur complement
159 
160 .vb
161    S = A11 - A10 inv(A00) A01
162 .ve
163 
164    PCLSC currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
165    inv(S) is given by
166 
167 .vb
168    inv(A10 A01) A10 A00 A01 inv(A10 A01)
169 .ve
170 
171    The product A10 A01 can be computed for you, but you can provide it (this is
172    usually more efficient anyway).  In the case of incompressible flow, A10 A10 is a Laplacian, call it L.  The current
173    interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.
174 
175    If you had called KSPSetOperators(ksp,S,Sp,flg), S should have type MATSCHURCOMPLEMENT and Sp can be any type you
176    like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
177    For example, you might have setup code like this
178 
179 .vb
180    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
181    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
182 .ve
183 
184    And then your Jacobian assembly would look like
185 
186 .vb
187    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
188    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
189    if (L) { assembly L }
190    if (Lp) { assemble Lp }
191 .ve
192 
193    With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L
194 
195 .vb
196    -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
197 .ve
198 
199    Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.
200 
201    Concepts: physics based preconditioners, block preconditioners
202 
203 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT,
204            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition(),
205            MatCreateSchurComplement()
206 M*/
207 
208 EXTERN_C_BEGIN
209 #undef __FUNCT__
210 #define __FUNCT__ "PCCreate_LSC"
211 PetscErrorCode  PCCreate_LSC(PC pc)
212 {
213   PC_LSC         *lsc;
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   ierr      = PetscNewLog(pc,PC_LSC,&lsc);CHKERRQ(ierr);
218   pc->data  = (void*)lsc;
219 
220   pc->ops->apply               = PCApply_LSC;
221   pc->ops->applytranspose      = 0;
222   pc->ops->setup               = PCSetUp_LSC;
223   pc->ops->destroy             = PCDestroy_LSC;
224   pc->ops->setfromoptions      = PCSetFromOptions_LSC;
225   pc->ops->view                = PCView_LSC;
226   pc->ops->applyrichardson     = 0;
227   PetscFunctionReturn(0);
228 }
229 EXTERN_C_END
230