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