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