xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 18be62a5feccf172f7bc80c15c4be8f6d6443e8b)
1 #define PETSCKSP_DLL
2 
3 /*
4       Defines a preconditioner that can consist of a collection of PCs
5 */
6 #include "src/ksp/pc/pcimpl.h"   /*I "petscpc.h" I*/
7 #include "petscksp.h"            /*I "petscksp.h" I*/
8 
9 typedef struct _PC_CompositeLink *PC_CompositeLink;
10 struct _PC_CompositeLink {
11   PC               pc;
12   PC_CompositeLink next;
13 };
14 
15 typedef struct {
16   PC_CompositeLink head;
17   PCCompositeType  type;
18   Vec              work1;
19   Vec              work2;
20   PetscScalar      alpha;
21   PetscTruth       use_true_matrix;
22 } PC_Composite;
23 
24 #undef __FUNCT__
25 #define __FUNCT__ "PCApply_Composite_Multiplicative"
26 static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
27 {
28   PetscErrorCode   ierr;
29   PC_Composite     *jac = (PC_Composite*)pc->data;
30   PC_CompositeLink next = jac->head;
31   PetscScalar      one = 1.0,mone = -1.0;
32   Mat              mat = pc->pmat;
33 
34   PetscFunctionBegin;
35   if (!next) {
36     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
37   }
38   if (next->next && !jac->work2) { /* allocate second work vector */
39     ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr);
40   }
41   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
42   if (jac->use_true_matrix) mat = pc->mat;
43   while (next->next) {
44     next = next->next;
45     ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr);
46     ierr = VecWAXPY(jac->work2,mone,jac->work1,x);CHKERRQ(ierr);
47     ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr);
48     ierr = VecAXPY(y,one,jac->work1);CHKERRQ(ierr);
49   }
50 
51   PetscFunctionReturn(0);
52 }
53 
54 /*
55     This is very special for a matrix of the form alpha I + R + S
56 where first preconditioner is built from alpha I + S and second from
57 alpha I + R
58 */
59 #undef __FUNCT__
60 #define __FUNCT__ "PCApply_Composite_Special"
61 static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
62 {
63   PetscErrorCode   ierr;
64   PC_Composite     *jac = (PC_Composite*)pc->data;
65   PC_CompositeLink next = jac->head;
66 
67   PetscFunctionBegin;
68   if (!next) {
69     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
70   }
71   if (!next->next || next->next->next) {
72     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
73   }
74 
75   ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
76   ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "PCApply_Composite_Additive"
82 static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
83 {
84   PetscErrorCode   ierr;
85   PC_Composite     *jac = (PC_Composite*)pc->data;
86   PC_CompositeLink next = jac->head;
87   PetscScalar      one = 1.0;
88 
89   PetscFunctionBegin;
90   if (!next) {
91     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
92   }
93   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
94   while (next->next) {
95     next = next->next;
96     ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
97     ierr = VecAXPY(y,one,jac->work1);CHKERRQ(ierr);
98   }
99   PetscFunctionReturn(0);
100 }
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "PCSetUp_Composite"
104 static PetscErrorCode PCSetUp_Composite(PC pc)
105 {
106   PetscErrorCode   ierr;
107   PC_Composite     *jac = (PC_Composite*)pc->data;
108   PC_CompositeLink next = jac->head;
109 
110   PetscFunctionBegin;
111   if (!jac->work1) {
112    ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr);
113   }
114   while (next) {
115     ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
116     next = next->next;
117   }
118 
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "PCDestroy_Composite"
124 static PetscErrorCode PCDestroy_Composite(PC pc)
125 {
126   PC_Composite     *jac = (PC_Composite*)pc->data;
127   PetscErrorCode   ierr;
128   PC_CompositeLink next = jac->head;
129 
130   PetscFunctionBegin;
131   while (next) {
132     ierr = PCDestroy(next->pc);CHKERRQ(ierr);
133     next = next->next;
134   }
135 
136   if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);}
137   if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);}
138   ierr = PetscFree(jac);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "PCSetFromOptions_Composite"
144 static PetscErrorCode PCSetFromOptions_Composite(PC pc)
145 {
146   PC_Composite     *jac = (PC_Composite*)pc->data;
147   PetscErrorCode   ierr;
148   PetscInt         nmax = 8,i;
149   PC_CompositeLink next;
150   char             *pcs[8];
151   PetscTruth       flg;
152 
153   PetscFunctionBegin;
154   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
155   ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
156     ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr);
157     if (flg) {
158       ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr);
159     }
160     ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr);
161     if (flg) {
162       for (i=0; i<nmax; i++) {
163         ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr);
164       }
165     }
166   ierr = PetscOptionsTail();CHKERRQ(ierr);
167 
168   next = jac->head;
169   while (next) {
170     ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr);
171     next = next->next;
172   }
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "PCView_Composite"
178 static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
179 {
180   PC_Composite     *jac = (PC_Composite*)pc->data;
181   PetscErrorCode   ierr;
182   PC_CompositeLink next = jac->head;
183   PetscTruth       iascii;
184 
185   PetscFunctionBegin;
186   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
187   if (iascii) {
188     ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr);
189     ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr);
190     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
191   } else {
192     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
193   }
194   if (iascii) {
195     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
196   }
197   while (next) {
198     ierr = PCView(next->pc,viewer);CHKERRQ(ierr);
199     next = next->next;
200   }
201   if (iascii) {
202     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
203     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
204   }
205   PetscFunctionReturn(0);
206 }
207 
208 /* ------------------------------------------------------------------------------*/
209 
210 EXTERN_C_BEGIN
211 #undef __FUNCT__
212 #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite"
213 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
214 {
215   PC_Composite *jac = (PC_Composite*)pc->data;
216   PetscFunctionBegin;
217   jac->alpha = alpha;
218   PetscFunctionReturn(0);
219 }
220 EXTERN_C_END
221 
222 EXTERN_C_BEGIN
223 #undef __FUNCT__
224 #define __FUNCT__ "PCCompositeSetType_Composite"
225 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type)
226 {
227   PetscFunctionBegin;
228   if (type == PC_COMPOSITE_ADDITIVE) {
229     pc->ops->apply = PCApply_Composite_Additive;
230   } else if (type ==  PC_COMPOSITE_MULTIPLICATIVE) {
231     pc->ops->apply = PCApply_Composite_Multiplicative;
232   } else if (type ==  PC_COMPOSITE_SPECIAL) {
233     pc->ops->apply = PCApply_Composite_Special;
234   } else {
235     SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
236   }
237   PetscFunctionReturn(0);
238 }
239 EXTERN_C_END
240 
241 EXTERN_C_BEGIN
242 #undef __FUNCT__
243 #define __FUNCT__ "PCCompositeAddPC_Composite"
244 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC_Composite(PC pc,PCType type)
245 {
246   PC_Composite     *jac;
247   PC_CompositeLink next,ilink;
248   PetscErrorCode   ierr;
249   PetscInt         cnt = 0;
250   const char       *prefix;
251   char             newprefix[8];
252 
253   PetscFunctionBegin;
254   ierr       = PetscNew(struct _PC_CompositeLink,&ilink);CHKERRQ(ierr);
255   ilink->next = 0;
256   ierr = PCCreate(pc->comm,&ilink->pc);CHKERRQ(ierr);
257 
258   jac  = (PC_Composite*)pc->data;
259   next = jac->head;
260   if (!next) {
261     jac->head = ilink;
262   } else {
263     cnt++;
264     while (next->next) {
265       next = next->next;
266       cnt++;
267     }
268     next->next = ilink;
269   }
270   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
271   ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr);
272   sprintf(newprefix,"sub_%d_",(int)cnt);
273   ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr);
274   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
275   ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr);
276 
277   PetscFunctionReturn(0);
278 }
279 EXTERN_C_END
280 
281 EXTERN_C_BEGIN
282 #undef __FUNCT__
283 #define __FUNCT__ "PCCompositeGetPC_Composite"
284 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
285 {
286   PC_Composite     *jac;
287   PC_CompositeLink next;
288   PetscInt         i;
289 
290   PetscFunctionBegin;
291   jac  = (PC_Composite*)pc->data;
292   next = jac->head;
293   for (i=0; i<n; i++) {
294     if (!next->next) {
295       SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
296     }
297     next = next->next;
298   }
299   *subpc = next->pc;
300   PetscFunctionReturn(0);
301 }
302 EXTERN_C_END
303 
304 EXTERN_C_BEGIN
305 #undef __FUNCT__
306 #define __FUNCT__ "PCCompositeSetUseTrue_Composite"
307 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc)
308 {
309   PC_Composite   *jac;
310 
311   PetscFunctionBegin;
312   jac                  = (PC_Composite*)pc->data;
313   jac->use_true_matrix = PETSC_TRUE;
314   PetscFunctionReturn(0);
315 }
316 EXTERN_C_END
317 
318 /* -------------------------------------------------------------------------------- */
319 #undef __FUNCT__
320 #define __FUNCT__ "PCCompositeSetType"
321 /*@C
322    PCCompositeSetType - Sets the type of composite preconditioner.
323 
324    Collective on PC
325 
326    Input Parameter:
327 .  pc - the preconditioner context
328 .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
329 
330    Options Database Key:
331 .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
332 
333    Level: Developer
334 
335 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
336 @*/
337 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type)
338 {
339   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
340 
341   PetscFunctionBegin;
342   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
343   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
344   if (f) {
345     ierr = (*f)(pc,type);CHKERRQ(ierr);
346   }
347   PetscFunctionReturn(0);
348 }
349 
350 #undef __FUNCT__
351 #define __FUNCT__ "PCCompositeSpecialSetAlpha"
352 /*@C
353    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
354      for alphaI + R + S
355 
356    Collective on PC
357 
358    Input Parameter:
359 +  pc - the preconditioner context
360 -  alpha - scale on identity
361 
362    Level: Developer
363 
364 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
365 @*/
366 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
367 {
368   PetscErrorCode ierr,(*f)(PC,PetscScalar);
369 
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
372   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr);
373   if (f) {
374     ierr = (*f)(pc,alpha);CHKERRQ(ierr);
375   }
376   PetscFunctionReturn(0);
377 }
378 
379 #undef __FUNCT__
380 #define __FUNCT__ "PCCompositeAddPC"
381 /*@C
382    PCCompositeAddPC - Adds another PC to the composite PC.
383 
384    Collective on PC
385 
386    Input Parameters:
387 .  pc - the preconditioner context
388 .  type - the type of the new preconditioner
389 
390    Level: Developer
391 
392 .keywords: PC, composite preconditioner, add
393 @*/
394 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type)
395 {
396   PetscErrorCode ierr,(*f)(PC,PCType);
397 
398   PetscFunctionBegin;
399   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
400   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr);
401   if (f) {
402     ierr = (*f)(pc,type);CHKERRQ(ierr);
403   }
404   PetscFunctionReturn(0);
405 }
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "PCCompositeGetPC"
409 /*@C
410    PCCompositeGetPC - Gets one of the PC objects in the composite PC.
411 
412    Not Collective
413 
414    Input Parameter:
415 .  pc - the preconditioner context
416 .  n - the number of the pc requested
417 
418    Output Parameters:
419 .  subpc - the PC requested
420 
421    Level: Developer
422 
423 .keywords: PC, get, composite preconditioner, sub preconditioner
424 
425 .seealso: PCCompositeAddPC()
426 @*/
427 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
428 {
429   PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
433   PetscValidPointer(subpc,3);
434   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
435   if (f) {
436     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
437   } else {
438     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
439   }
440   PetscFunctionReturn(0);
441 }
442 
443 #undef __FUNCT__
444 #define __FUNCT__ "PCCompositeSetUseTrue"
445 /*@
446    PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
447                       the matrix used to define the preconditioner) is used to compute
448                       the residual when the multiplicative scheme is used.
449 
450    Collective on PC
451 
452    Input Parameters:
453 .  pc - the preconditioner context
454 
455    Options Database Key:
456 .  -pc_composite_true - Activates PCCompositeSetUseTrue()
457 
458    Note:
459    For the common case in which the preconditioning and linear
460    system matrices are identical, this routine is unnecessary.
461 
462    Level: Developer
463 
464 .keywords: PC, composite preconditioner, set, true, flag
465 
466 .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
467 @*/
468 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc)
469 {
470   PetscErrorCode ierr,(*f)(PC);
471 
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
474   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr);
475   if (f) {
476     ierr = (*f)(pc);CHKERRQ(ierr);
477   }
478   PetscFunctionReturn(0);
479 }
480 
481 /* -------------------------------------------------------------------------------------------*/
482 
483 /*MC
484      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
485 
486    Options Database Keys:
487 .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
488 .  -pc_composite_true - Activates PCCompositeSetUseTrue()
489 
490    Level: intermediate
491 
492    Concepts: composing solvers
493 
494    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
495           inner PCs to be PCKSP.
496           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
497           the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
498 
499 
500 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
501            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
502            PCCompositeGetPC(), PCCompositeSetUseTrue()
503 
504 M*/
505 
506 EXTERN_C_BEGIN
507 #undef __FUNCT__
508 #define __FUNCT__ "PCCreate_Composite"
509 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc)
510 {
511   PetscErrorCode ierr;
512   PC_Composite   *jac;
513 
514   PetscFunctionBegin;
515   ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr);
516   ierr = PetscLogObjectMemory(pc,sizeof(PC_Composite));CHKERRQ(ierr);
517   pc->ops->apply              = PCApply_Composite_Additive;
518   pc->ops->setup              = PCSetUp_Composite;
519   pc->ops->destroy            = PCDestroy_Composite;
520   pc->ops->setfromoptions     = PCSetFromOptions_Composite;
521   pc->ops->view               = PCView_Composite;
522   pc->ops->applyrichardson    = 0;
523 
524   pc->data               = (void*)jac;
525   jac->type              = PC_COMPOSITE_ADDITIVE;
526   jac->work1             = 0;
527   jac->work2             = 0;
528   jac->head              = 0;
529 
530   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
531                     PCCompositeSetType_Composite);CHKERRQ(ierr);
532   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
533                     PCCompositeAddPC_Composite);CHKERRQ(ierr);
534   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
535                     PCCompositeGetPC_Composite);CHKERRQ(ierr);
536   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
537                     PCCompositeSetUseTrue_Composite);CHKERRQ(ierr);
538   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
539                     PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr);
540 
541   PetscFunctionReturn(0);
542 }
543 EXTERN_C_END
544 
545