xref: /petsc/src/ksp/pc/impls/galerkin/galerkin.c (revision 90e035378b7ae2ceaeb17618b77f7bb7d863b669)
1 
2 /*
3       Defines a preconditioner defined by R^T S R
4 */
5 #include <petsc/private/pcimpl.h>
6 #include <petscksp.h>         /*I "petscksp.h" I*/
7 
8 typedef struct {
9   KSP            ksp;
10   Mat            R,P;
11   Vec            b,x;
12   PetscErrorCode (*computeasub)(PC,Mat,Mat,Mat*,void*);
13   void           *computeasub_ctx;
14 } PC_Galerkin;
15 
16 static PetscErrorCode PCApply_Galerkin(PC pc,Vec x,Vec y)
17 {
18   PetscErrorCode ierr;
19   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
20 
21   PetscFunctionBegin;
22   if (jac->R) {
23     ierr = MatRestrict(jac->R,x,jac->b);CHKERRQ(ierr);
24   } else {
25     ierr = MatRestrict(jac->P,x,jac->b);CHKERRQ(ierr);
26   }
27   ierr = KSPSolve(jac->ksp,jac->b,jac->x);CHKERRQ(ierr);
28   ierr = KSPCheckSolve(jac->ksp,pc,jac->x);CHKERRQ(ierr);
29   if (jac->P) {
30     ierr = MatInterpolate(jac->P,jac->x,y);CHKERRQ(ierr);
31   } else {
32     ierr = MatInterpolate(jac->R,jac->x,y);CHKERRQ(ierr);
33   }
34   PetscFunctionReturn(0);
35 }
36 
37 static PetscErrorCode PCSetUp_Galerkin(PC pc)
38 {
39   PetscErrorCode ierr;
40   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
41   PetscBool      a;
42   Vec            *xx,*yy;
43 
44   PetscFunctionBegin;
45   if (jac->computeasub) {
46     Mat Ap;
47     if (!pc->setupcalled) {
48       ierr = (*jac->computeasub)(pc,pc->pmat,NULL,&Ap,jac->computeasub_ctx);CHKERRQ(ierr);
49       ierr = KSPSetOperators(jac->ksp,Ap,Ap);CHKERRQ(ierr);
50       ierr = MatDestroy(&Ap);CHKERRQ(ierr);
51     } else {
52       ierr = KSPGetOperators(jac->ksp,NULL,&Ap);CHKERRQ(ierr);
53       ierr = (*jac->computeasub)(pc,pc->pmat,Ap,NULL,jac->computeasub_ctx);CHKERRQ(ierr);
54     }
55   }
56 
57   if (!jac->x) {
58     ierr = KSPGetOperatorsSet(jac->ksp,&a,NULL);CHKERRQ(ierr);
59     if (!a) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set operator of PCGALERKIN KSP with PCGalerkinGetKSP()/KSPSetOperators()");
60     ierr   = KSPCreateVecs(jac->ksp,1,&xx,1,&yy);CHKERRQ(ierr);
61     jac->x = *xx;
62     jac->b = *yy;
63     ierr   = PetscFree(xx);CHKERRQ(ierr);
64     ierr   = PetscFree(yy);CHKERRQ(ierr);
65   }
66   if (!jac->R && !jac->P) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or interpolation of PCGALERKIN with PCGalerkinSetRestriction()/Interpolation()");
67   /* should check here that sizes of R/P match size of a */
68 
69   PetscFunctionReturn(0);
70 }
71 
72 static PetscErrorCode PCReset_Galerkin(PC pc)
73 {
74   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   ierr = MatDestroy(&jac->R);CHKERRQ(ierr);
79   ierr = MatDestroy(&jac->P);CHKERRQ(ierr);
80   ierr = VecDestroy(&jac->x);CHKERRQ(ierr);
81   ierr = VecDestroy(&jac->b);CHKERRQ(ierr);
82   ierr = KSPReset(jac->ksp);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 static PetscErrorCode PCDestroy_Galerkin(PC pc)
87 {
88   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
89   PetscErrorCode ierr;
90 
91   PetscFunctionBegin;
92   ierr = PCReset_Galerkin(pc);CHKERRQ(ierr);
93   ierr = KSPDestroy(&jac->ksp);CHKERRQ(ierr);
94   ierr = PetscFree(pc->data);CHKERRQ(ierr);
95   PetscFunctionReturn(0);
96 }
97 
98 static PetscErrorCode PCView_Galerkin(PC pc,PetscViewer viewer)
99 {
100   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
101   PetscErrorCode ierr;
102   PetscBool      iascii;
103 
104   PetscFunctionBegin;
105   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
106   if (iascii) {
107     ierr = PetscViewerASCIIPrintf(viewer,"  KSP on Galerkin follow\n");CHKERRQ(ierr);
108     ierr = PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n");CHKERRQ(ierr);
109   }
110   ierr = KSPView(jac->ksp,viewer);CHKERRQ(ierr);
111   PetscFunctionReturn(0);
112 }
113 
114 static PetscErrorCode  PCGalerkinGetKSP_Galerkin(PC pc,KSP *ksp)
115 {
116   PC_Galerkin *jac = (PC_Galerkin*)pc->data;
117 
118   PetscFunctionBegin;
119   *ksp = jac->ksp;
120   PetscFunctionReturn(0);
121 }
122 
123 static PetscErrorCode  PCGalerkinSetRestriction_Galerkin(PC pc,Mat R)
124 {
125   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129   ierr   = PetscObjectReference((PetscObject)R);CHKERRQ(ierr);
130   ierr   = MatDestroy(&jac->R);CHKERRQ(ierr);
131   jac->R = R;
132   PetscFunctionReturn(0);
133 }
134 
135 static PetscErrorCode  PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P)
136 {
137   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   ierr   = PetscObjectReference((PetscObject)P);CHKERRQ(ierr);
142   ierr   = MatDestroy(&jac->P);CHKERRQ(ierr);
143   jac->P = P;
144   PetscFunctionReturn(0);
145 }
146 
147 static PetscErrorCode  PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx)
148 {
149   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
150 
151   PetscFunctionBegin;
152   jac->computeasub     = computeAsub;
153   jac->computeasub_ctx = ctx;
154   PetscFunctionReturn(0);
155 }
156 
157 /* -------------------------------------------------------------------------------- */
158 /*@
159    PCGalerkinSetRestriction - Sets the restriction operator for the "Galerkin-type" preconditioner
160 
161    Logically Collective on PC
162 
163    Input Parameter:
164 +  pc - the preconditioner context
165 -  R - the restriction operator
166 
167    Notes:
168     Either this or PCGalerkinSetInterpolation() or both must be called
169 
170    Level: Intermediate
171 
172 .keywords: PC, set, Galerkin preconditioner
173 
174 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
175            PCGalerkinSetInterpolation(), PCGalerkinGetKSP()
176 
177 @*/
178 PetscErrorCode  PCGalerkinSetRestriction(PC pc,Mat R)
179 {
180   PetscErrorCode ierr;
181 
182   PetscFunctionBegin;
183   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
184   ierr = PetscTryMethod(pc,"PCGalerkinSetRestriction_C",(PC,Mat),(pc,R));CHKERRQ(ierr);
185   PetscFunctionReturn(0);
186 }
187 
188 /*@
189    PCGalerkinSetInterpolation - Sets the interpolation operator for the "Galerkin-type" preconditioner
190 
191    Logically Collective on PC
192 
193    Input Parameter:
194 +  pc - the preconditioner context
195 -  R - the interpolation operator
196 
197    Notes:
198     Either this or PCGalerkinSetRestriction() or both must be called
199 
200    Level: Intermediate
201 
202 .keywords: PC, set, Galerkin preconditioner
203 
204 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
205            PCGalerkinSetRestriction(), PCGalerkinGetKSP()
206 
207 @*/
208 PetscErrorCode  PCGalerkinSetInterpolation(PC pc,Mat P)
209 {
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
214   ierr = PetscTryMethod(pc,"PCGalerkinSetInterpolation_C",(PC,Mat),(pc,P));CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 /*@
219    PCGalerkinSetComputeSubmatrix - Provide a routine that will be called to compute the Galerkin submatrix
220 
221    Logically Collective
222 
223    Input Parameter:
224 +  pc - the preconditioner context
225 .  computeAsub - routine that computes the submatrix from the global matrix
226 -  ctx - context used by the routine, or NULL
227 
228    Calling sequence of computeAsub:
229 $    computeAsub(PC pc,Mat A, Mat Ap, Mat *cAP,void *ctx);
230 
231 +  PC - the Galerkin PC
232 .  A - the matrix in the Galerkin PC
233 .  Ap - the computed submatrix from any previous computation, if NULL it has not previously been computed
234 .  cAp - the submatrix computed by this routine
235 -  ctx - optional user-defined function context
236 
237    Level: Intermediate
238 
239    Notes:
240     Instead of providing this routine you can call PCGalerkinGetKSP() and then KSPSetOperators() to provide the submatrix,
241           but that will not work for multiple KSPSolves with different matrices unless you call it for each solve.
242 
243           This routine is called each time the outter matrix is changed. In the first call the Ap argument is NULL and the routine should create the
244           matrix and computes its values in cAp. On each subsequent call the routine should up the Ap matrix.
245 
246    Developer Notes:
247     If the user does not call this routine nor call PCGalerkinGetKSP() and KSPSetOperators() then PCGalerkin could
248                     could automatically compute the submatrix via calls to MatGalerkin() or MatRARt()
249 
250 .keywords: PC, get, Galerkin preconditioner, sub preconditioner
251 
252 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
253            PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP()
254 
255 @*/
256 PetscErrorCode  PCGalerkinSetComputeSubmatrix(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx)
257 {
258   PetscErrorCode ierr;
259 
260   PetscFunctionBegin;
261   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
262   ierr = PetscTryMethod(pc,"PCGalerkinSetComputeSubmatrix_C",(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*),(pc,computeAsub,ctx));CHKERRQ(ierr);
263   PetscFunctionReturn(0);
264 }
265 
266 /*@
267    PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC.
268 
269    Not Collective
270 
271    Input Parameter:
272 .  pc - the preconditioner context
273 
274    Output Parameters:
275 .  ksp - the KSP object
276 
277    Level: Intermediate
278 
279    Notes:
280     Once you have called this routine you can call KSPSetOperators() on the resulting ksp to provide the operator for the Galerkin problem,
281           an alternative is to use PCGalerkinSetComputeSubmatrix() to provide a routine that computes the submatrix as needed.
282 
283 .keywords: PC, get, Galerkin preconditioner, sub preconditioner
284 
285 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
286            PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinSetComputeSubmatrix()
287 
288 @*/
289 PetscErrorCode  PCGalerkinGetKSP(PC pc,KSP *ksp)
290 {
291   PetscErrorCode ierr;
292 
293   PetscFunctionBegin;
294   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
295   PetscValidPointer(ksp,2);
296   ierr = PetscUseMethod(pc,"PCGalerkinGetKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 static PetscErrorCode PCSetFromOptions_Galerkin(PetscOptionItems *PetscOptionsObject,PC pc)
301 {
302   PC_Galerkin    *jac = (PC_Galerkin*)pc->data;
303   PetscErrorCode ierr;
304   const char     *prefix;
305   PetscBool      flg;
306 
307   PetscFunctionBegin;
308   ierr = KSPGetOptionsPrefix(jac->ksp,&prefix);CHKERRQ(ierr);
309   ierr = PetscStrendswith(prefix,"galerkin_",&flg);CHKERRQ(ierr);
310   if (!flg) {
311     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
312     ierr = KSPSetOptionsPrefix(jac->ksp,prefix);CHKERRQ(ierr);
313     ierr = KSPAppendOptionsPrefix(jac->ksp,"galerkin_");CHKERRQ(ierr);
314   }
315 
316   ierr = PetscOptionsHead(PetscOptionsObject,"Galerkin options");CHKERRQ(ierr);
317   if (jac->ksp) {
318     ierr = KSPSetFromOptions(jac->ksp);CHKERRQ(ierr);
319   }
320   ierr = PetscOptionsTail();CHKERRQ(ierr);
321   PetscFunctionReturn(0);
322 }
323 
324 /* -------------------------------------------------------------------------------------------*/
325 
326 /*MC
327      PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T)
328 
329 $   Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by
330 $   PCGalerkinGetKSP(pc,&ksp); KSPSetOperators(ksp,A,....)
331 
332    Level: intermediate
333 
334    Developer Note: If KSPSetOperators() has not been called on the inner KSP then PCGALERKIN could use MatRARt() or MatPtAP() to compute
335                    the operators automatically.
336                    Should there be a prefix for the inner KSP.
337                    There is no KSPSetFromOptions_Galerkin() that calls KSPSetFromOptions() on the inner KSP
338 
339 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
340            PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP()
341 
342 M*/
343 
344 PETSC_EXTERN PetscErrorCode PCCreate_Galerkin(PC pc)
345 {
346   PetscErrorCode ierr;
347   PC_Galerkin    *jac;
348 
349   PetscFunctionBegin;
350   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
351 
352   pc->ops->apply           = PCApply_Galerkin;
353   pc->ops->setup           = PCSetUp_Galerkin;
354   pc->ops->reset           = PCReset_Galerkin;
355   pc->ops->destroy         = PCDestroy_Galerkin;
356   pc->ops->view            = PCView_Galerkin;
357   pc->ops->setfromoptions  = PCSetFromOptions_Galerkin;
358   pc->ops->applyrichardson = NULL;
359 
360   ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->ksp);CHKERRQ(ierr);
361   ierr = KSPSetErrorIfNotConverged(jac->ksp,pc->erroriffailure);CHKERRQ(ierr);
362   ierr = PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
363 
364   pc->data = (void*)jac;
365 
366   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetRestriction_C",PCGalerkinSetRestriction_Galerkin);CHKERRQ(ierr);
367   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetInterpolation_C",PCGalerkinSetInterpolation_Galerkin);CHKERRQ(ierr);
368   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinGetKSP_C",PCGalerkinGetKSP_Galerkin);CHKERRQ(ierr);
369   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetComputeSubmatrix_C",PCGalerkinSetComputeSubmatrix_Galerkin);CHKERRQ(ierr);
370   PetscFunctionReturn(0);
371 }
372 
373