xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 97c4aaa0041b17e0788ef5380763993e79c96a23)
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   PetscBool        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   PetscBool        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   PetscBool        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   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->pc);CHKERRQ(ierr);
264 
265   jac  = (PC_Composite*)pc->data;
266   next = jac->head;
267   if (!next) {
268     jac->head       = ilink;
269     ilink->previous = PETSC_NULL;
270   } else {
271     cnt++;
272     while (next->next) {
273       next = next->next;
274       cnt++;
275     }
276     next->next      = ilink;
277     ilink->previous = next;
278   }
279   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
280   ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr);
281   sprintf(newprefix,"sub_%d_",(int)cnt);
282   ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr);
283   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
284   ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 EXTERN_C_END
288 
289 EXTERN_C_BEGIN
290 #undef __FUNCT__
291 #define __FUNCT__ "PCCompositeGetPC_Composite"
292 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
293 {
294   PC_Composite     *jac;
295   PC_CompositeLink next;
296   PetscInt         i;
297 
298   PetscFunctionBegin;
299   jac  = (PC_Composite*)pc->data;
300   next = jac->head;
301   for (i=0; i<n; i++) {
302     if (!next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
303     next = next->next;
304   }
305   *subpc = next->pc;
306   PetscFunctionReturn(0);
307 }
308 EXTERN_C_END
309 
310 EXTERN_C_BEGIN
311 #undef __FUNCT__
312 #define __FUNCT__ "PCCompositeSetUseTrue_Composite"
313 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc)
314 {
315   PC_Composite   *jac;
316 
317   PetscFunctionBegin;
318   jac                  = (PC_Composite*)pc->data;
319   jac->use_true_matrix = PETSC_TRUE;
320   PetscFunctionReturn(0);
321 }
322 EXTERN_C_END
323 
324 /* -------------------------------------------------------------------------------- */
325 #undef __FUNCT__
326 #define __FUNCT__ "PCCompositeSetType"
327 /*@
328    PCCompositeSetType - Sets the type of composite preconditioner.
329 
330    Logically Collective on PC
331 
332    Input Parameter:
333 +  pc - the preconditioner context
334 -  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
335 
336    Options Database Key:
337 .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
338 
339    Level: Developer
340 
341 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
342 @*/
343 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type)
344 {
345   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
346 
347   PetscFunctionBegin;
348   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
349   PetscValidLogicalCollectiveEnum(pc,type,2);
350 
351   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
352   if (f) {
353     ierr = (*f)(pc,type);CHKERRQ(ierr);
354   }
355   PetscFunctionReturn(0);
356 }
357 
358 #undef __FUNCT__
359 #define __FUNCT__ "PCCompositeSpecialSetAlpha"
360 /*@
361    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
362      for alphaI + R + S
363 
364    Logically Collective on PC
365 
366    Input Parameter:
367 +  pc - the preconditioner context
368 -  alpha - scale on identity
369 
370    Level: Developer
371 
372 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
373 @*/
374 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
375 {
376   PetscErrorCode ierr,(*f)(PC,PetscScalar);
377 
378   PetscFunctionBegin;
379   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
380   PetscValidLogicalCollectiveScalar(pc,alpha,2);
381 
382   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr);
383   if (f) {
384     ierr = (*f)(pc,alpha);CHKERRQ(ierr);
385   }
386   PetscFunctionReturn(0);
387 }
388 
389 #undef __FUNCT__
390 #define __FUNCT__ "PCCompositeAddPC"
391 /*@C
392    PCCompositeAddPC - Adds another PC to the composite PC.
393 
394    Collective on PC
395 
396    Input Parameters:
397 +  pc - the preconditioner context
398 -  type - the type of the new preconditioner
399 
400    Level: Developer
401 
402 .keywords: PC, composite preconditioner, add
403 @*/
404 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type)
405 {
406   PetscErrorCode ierr,(*f)(PC,PCType);
407 
408   PetscFunctionBegin;
409   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
410   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr);
411   if (f) {
412     ierr = (*f)(pc,type);CHKERRQ(ierr);
413   }
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "PCCompositeGetPC"
419 /*@
420    PCCompositeGetPC - Gets one of the PC objects in the composite PC.
421 
422    Not Collective
423 
424    Input Parameter:
425 +  pc - the preconditioner context
426 -  n - the number of the pc requested
427 
428    Output Parameters:
429 .  subpc - the PC requested
430 
431    Level: Developer
432 
433 .keywords: PC, get, composite preconditioner, sub preconditioner
434 
435 .seealso: PCCompositeAddPC()
436 @*/
437 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
438 {
439   PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
440 
441   PetscFunctionBegin;
442   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
443   PetscValidPointer(subpc,3);
444   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
445   if (f) {
446     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
447   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "PCCompositeSetUseTrue"
453 /*@
454    PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
455                       the matrix used to define the preconditioner) is used to compute
456                       the residual when the multiplicative scheme is used.
457 
458    Logically Collective on PC
459 
460    Input Parameters:
461 .  pc - the preconditioner context
462 
463    Options Database Key:
464 .  -pc_composite_true - Activates PCCompositeSetUseTrue()
465 
466    Note:
467    For the common case in which the preconditioning and linear
468    system matrices are identical, this routine is unnecessary.
469 
470    Level: Developer
471 
472 .keywords: PC, composite preconditioner, set, true, flag
473 
474 .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
475 @*/
476 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc)
477 {
478   PetscErrorCode ierr,(*f)(PC);
479 
480   PetscFunctionBegin;
481   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
482   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr);
483   if (f) {
484     ierr = (*f)(pc);CHKERRQ(ierr);
485   }
486   PetscFunctionReturn(0);
487 }
488 
489 /* -------------------------------------------------------------------------------------------*/
490 
491 /*MC
492      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
493 
494    Options Database Keys:
495 +  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
496 .  -pc_composite_true - Activates PCCompositeSetUseTrue()
497 -  -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose
498 
499    Level: intermediate
500 
501    Concepts: composing solvers
502 
503    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
504           inner PCs to be PCKSP.
505           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
506           the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
507 
508 
509 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
510            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
511            PCCompositeGetPC(), PCCompositeSetUseTrue()
512 
513 M*/
514 
515 EXTERN_C_BEGIN
516 #undef __FUNCT__
517 #define __FUNCT__ "PCCreate_Composite"
518 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc)
519 {
520   PetscErrorCode ierr;
521   PC_Composite   *jac;
522 
523   PetscFunctionBegin;
524   ierr = PetscNewLog(pc,PC_Composite,&jac);CHKERRQ(ierr);
525   pc->ops->apply              = PCApply_Composite_Additive;
526   pc->ops->setup              = PCSetUp_Composite;
527   pc->ops->destroy            = PCDestroy_Composite;
528   pc->ops->setfromoptions     = PCSetFromOptions_Composite;
529   pc->ops->view               = PCView_Composite;
530   pc->ops->applyrichardson    = 0;
531 
532   pc->data               = (void*)jac;
533   jac->type              = PC_COMPOSITE_ADDITIVE;
534   jac->work1             = 0;
535   jac->work2             = 0;
536   jac->head              = 0;
537 
538   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
539                     PCCompositeSetType_Composite);CHKERRQ(ierr);
540   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
541                     PCCompositeAddPC_Composite);CHKERRQ(ierr);
542   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
543                     PCCompositeGetPC_Composite);CHKERRQ(ierr);
544   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
545                     PCCompositeSetUseTrue_Composite);CHKERRQ(ierr);
546   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
547                     PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr);
548 
549   PetscFunctionReturn(0);
550 }
551 EXTERN_C_END
552 
553