xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 5f3c5e7ab1713b2b36ec2007ece43899b4f0dcb3)
1 
2 /*
3       Defines a preconditioner that can consist of a collection of PCs
4 */
5 #include <petsc/private/pcimpl.h>
6 #include <petscksp.h> /*I "petscksp.h" I*/
7 
8 typedef struct _PC_CompositeLink *PC_CompositeLink;
9 struct _PC_CompositeLink {
10   PC               pc;
11   PC_CompositeLink next;
12   PC_CompositeLink previous;
13 };
14 
15 typedef struct {
16   PC_CompositeLink head;
17   PCCompositeType  type;
18   Vec              work1;
19   Vec              work2;
20   PetscScalar      alpha;
21 } PC_Composite;
22 
23 static PetscErrorCode PCApply_Composite_Multiplicative(PC pc, Vec x, Vec y)
24 {
25   PC_Composite    *jac  = (PC_Composite *)pc->data;
26   PC_CompositeLink next = jac->head;
27   Mat              mat  = pc->pmat;
28 
29   PetscFunctionBegin;
30 
31   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
32 
33   /* Set the reuse flag on children PCs */
34   while (next) {
35     PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner));
36     next = next->next;
37   }
38   next = jac->head;
39 
40   if (next->next && !jac->work2) { /* allocate second work vector */
41     PetscCall(VecDuplicate(jac->work1, &jac->work2));
42   }
43   if (pc->useAmat) mat = pc->mat;
44   PetscCall(PCApply(next->pc, x, y)); /* y <- B x */
45   while (next->next) {
46     next = next->next;
47     PetscCall(MatMult(mat, y, jac->work1));               /* work1 <- A y */
48     PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x)); /* work2 <- x - work1 */
49     PetscCall(PCApply(next->pc, jac->work2, jac->work1)); /* work1 <- C work2 */
50     PetscCall(VecAXPY(y, 1.0, jac->work1));               /* y <- y + work1 = B x + C (x - A B x) = (B + C (1 - A B)) x */
51   }
52   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
53     while (next->previous) {
54       next = next->previous;
55       PetscCall(MatMult(mat, y, jac->work1));
56       PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x));
57       PetscCall(PCApply(next->pc, jac->work2, jac->work1));
58       PetscCall(VecAXPY(y, 1.0, jac->work1));
59     }
60   }
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63 
64 static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc, Vec x, Vec y)
65 {
66   PC_Composite    *jac  = (PC_Composite *)pc->data;
67   PC_CompositeLink next = jac->head;
68   Mat              mat  = pc->pmat;
69 
70   PetscFunctionBegin;
71   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
72   if (next->next && !jac->work2) { /* allocate second work vector */
73     PetscCall(VecDuplicate(jac->work1, &jac->work2));
74   }
75   if (pc->useAmat) mat = pc->mat;
76   /* locate last PC */
77   while (next->next) next = next->next;
78   PetscCall(PCApplyTranspose(next->pc, x, y));
79   while (next->previous) {
80     next = next->previous;
81     PetscCall(MatMultTranspose(mat, y, jac->work1));
82     PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x));
83     PetscCall(PCApplyTranspose(next->pc, jac->work2, jac->work1));
84     PetscCall(VecAXPY(y, 1.0, jac->work1));
85   }
86   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
87     next = jac->head;
88     while (next->next) {
89       next = next->next;
90       PetscCall(MatMultTranspose(mat, y, jac->work1));
91       PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x));
92       PetscCall(PCApplyTranspose(next->pc, jac->work2, jac->work1));
93       PetscCall(VecAXPY(y, 1.0, jac->work1));
94     }
95   }
96   PetscFunctionReturn(PETSC_SUCCESS);
97 }
98 
99 /*
100     This is very special for a matrix of the form alpha I + R + S
101 where first preconditioner is built from alpha I + S and second from
102 alpha I + R
103 */
104 static PetscErrorCode PCApply_Composite_Special(PC pc, Vec x, Vec y)
105 {
106   PC_Composite    *jac  = (PC_Composite *)pc->data;
107   PC_CompositeLink next = jac->head;
108 
109   PetscFunctionBegin;
110   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
111   PetscCheck(next->next && !next->next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Special composite preconditioners requires exactly two PCs");
112 
113   /* Set the reuse flag on children PCs */
114   PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner));
115   PetscCall(PCSetReusePreconditioner(next->next->pc, pc->reusepreconditioner));
116 
117   PetscCall(PCApply(next->pc, x, jac->work1));
118   PetscCall(PCApply(next->next->pc, jac->work1, y));
119   PetscFunctionReturn(PETSC_SUCCESS);
120 }
121 
122 static PetscErrorCode PCApply_Composite_Additive(PC pc, Vec x, Vec y)
123 {
124   PC_Composite    *jac  = (PC_Composite *)pc->data;
125   PC_CompositeLink next = jac->head;
126 
127   PetscFunctionBegin;
128   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
129 
130   /* Set the reuse flag on children PCs */
131   while (next) {
132     PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner));
133     next = next->next;
134   }
135   next = jac->head;
136 
137   PetscCall(PCApply(next->pc, x, y));
138   while (next->next) {
139     next = next->next;
140     PetscCall(PCApply(next->pc, x, jac->work1));
141     PetscCall(VecAXPY(y, 1.0, jac->work1));
142   }
143   PetscFunctionReturn(PETSC_SUCCESS);
144 }
145 
146 static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc, Vec x, Vec y)
147 {
148   PC_Composite    *jac  = (PC_Composite *)pc->data;
149   PC_CompositeLink next = jac->head;
150 
151   PetscFunctionBegin;
152   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
153   PetscCall(PCApplyTranspose(next->pc, x, y));
154   while (next->next) {
155     next = next->next;
156     PetscCall(PCApplyTranspose(next->pc, x, jac->work1));
157     PetscCall(VecAXPY(y, 1.0, jac->work1));
158   }
159   PetscFunctionReturn(PETSC_SUCCESS);
160 }
161 
162 static PetscErrorCode PCSetUp_Composite(PC pc)
163 {
164   PC_Composite    *jac  = (PC_Composite *)pc->data;
165   PC_CompositeLink next = jac->head;
166   DM               dm;
167 
168   PetscFunctionBegin;
169   if (!jac->work1) PetscCall(MatCreateVecs(pc->pmat, &jac->work1, NULL));
170   PetscCall(PCGetDM(pc, &dm));
171   while (next) {
172     if (!next->pc->dm) PetscCall(PCSetDM(next->pc, dm));
173     if (!next->pc->mat) PetscCall(PCSetOperators(next->pc, pc->mat, pc->pmat));
174     next = next->next;
175   }
176   PetscFunctionReturn(PETSC_SUCCESS);
177 }
178 
179 static PetscErrorCode PCSetUpOnBlocks_Composite(PC pc)
180 {
181   PC_Composite    *jac  = (PC_Composite *)pc->data;
182   PC_CompositeLink next = jac->head;
183   PCFailedReason   reason;
184 
185   PetscFunctionBegin;
186   while (next) {
187     PetscCall(PCSetUp(next->pc));
188     PetscCall(PCGetFailedReasonRank(next->pc, &reason));
189     if (reason) pc->failedreason = reason;
190     next = next->next;
191   }
192   PetscFunctionReturn(PETSC_SUCCESS);
193 }
194 
195 static PetscErrorCode PCReset_Composite(PC pc)
196 {
197   PC_Composite    *jac  = (PC_Composite *)pc->data;
198   PC_CompositeLink next = jac->head;
199 
200   PetscFunctionBegin;
201   while (next) {
202     PetscCall(PCReset(next->pc));
203     next = next->next;
204   }
205   PetscCall(VecDestroy(&jac->work1));
206   PetscCall(VecDestroy(&jac->work2));
207   PetscFunctionReturn(PETSC_SUCCESS);
208 }
209 
210 static PetscErrorCode PCDestroy_Composite(PC pc)
211 {
212   PC_Composite    *jac  = (PC_Composite *)pc->data;
213   PC_CompositeLink next = jac->head, next_tmp;
214 
215   PetscFunctionBegin;
216   PetscCall(PCReset_Composite(pc));
217   while (next) {
218     PetscCall(PCDestroy(&next->pc));
219     next_tmp = next;
220     next     = next->next;
221     PetscCall(PetscFree(next_tmp));
222   }
223   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", NULL));
224   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", NULL));
225   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", NULL));
226   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", NULL));
227   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", NULL));
228   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", NULL));
229   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", NULL));
230   PetscCall(PetscFree(pc->data));
231   PetscFunctionReturn(PETSC_SUCCESS);
232 }
233 
234 static PetscErrorCode PCSetFromOptions_Composite(PC pc, PetscOptionItems *PetscOptionsObject)
235 {
236   PC_Composite    *jac  = (PC_Composite *)pc->data;
237   PetscInt         nmax = 8, i;
238   PC_CompositeLink next;
239   char            *pcs[8];
240   PetscBool        flg;
241 
242   PetscFunctionBegin;
243   PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options");
244   PetscCall(PetscOptionsEnum("-pc_composite_type", "Type of composition", "PCCompositeSetType", PCCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg));
245   if (flg) PetscCall(PCCompositeSetType(pc, jac->type));
246   PetscCall(PetscOptionsStringArray("-pc_composite_pcs", "List of composite solvers", "PCCompositeAddPCType", pcs, &nmax, &flg));
247   if (flg) {
248     for (i = 0; i < nmax; i++) {
249       PetscCall(PCCompositeAddPCType(pc, pcs[i]));
250       PetscCall(PetscFree(pcs[i])); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */
251     }
252   }
253   PetscOptionsHeadEnd();
254 
255   next = jac->head;
256   while (next) {
257     PetscCall(PCSetFromOptions(next->pc));
258     next = next->next;
259   }
260   PetscFunctionReturn(PETSC_SUCCESS);
261 }
262 
263 static PetscErrorCode PCView_Composite(PC pc, PetscViewer viewer)
264 {
265   PC_Composite    *jac  = (PC_Composite *)pc->data;
266   PC_CompositeLink next = jac->head;
267   PetscBool        iascii;
268 
269   PetscFunctionBegin;
270   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
271   if (iascii) {
272     PetscCall(PetscViewerASCIIPrintf(viewer, "Composite PC type - %s\n", PCCompositeTypes[jac->type]));
273     PetscCall(PetscViewerASCIIPrintf(viewer, "PCs on composite preconditioner follow\n"));
274     PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n"));
275   }
276   if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer));
277   while (next) {
278     PetscCall(PCView(next->pc, viewer));
279     next = next->next;
280   }
281   if (iascii) {
282     PetscCall(PetscViewerASCIIPopTab(viewer));
283     PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n"));
284   }
285   PetscFunctionReturn(PETSC_SUCCESS);
286 }
287 
288 static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc, PetscScalar alpha)
289 {
290   PC_Composite *jac = (PC_Composite *)pc->data;
291 
292   PetscFunctionBegin;
293   jac->alpha = alpha;
294   PetscFunctionReturn(PETSC_SUCCESS);
295 }
296 
297 static PetscErrorCode PCCompositeSetType_Composite(PC pc, PCCompositeType type)
298 {
299   PC_Composite *jac = (PC_Composite *)pc->data;
300 
301   PetscFunctionBegin;
302   if (type == PC_COMPOSITE_ADDITIVE) {
303     pc->ops->apply          = PCApply_Composite_Additive;
304     pc->ops->applytranspose = PCApplyTranspose_Composite_Additive;
305   } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
306     pc->ops->apply          = PCApply_Composite_Multiplicative;
307     pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative;
308   } else if (type == PC_COMPOSITE_SPECIAL) {
309     pc->ops->apply          = PCApply_Composite_Special;
310     pc->ops->applytranspose = NULL;
311   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Unknown composite preconditioner type");
312   jac->type = type;
313   PetscFunctionReturn(PETSC_SUCCESS);
314 }
315 
316 static PetscErrorCode PCCompositeGetType_Composite(PC pc, PCCompositeType *type)
317 {
318   PC_Composite *jac = (PC_Composite *)pc->data;
319 
320   PetscFunctionBegin;
321   *type = jac->type;
322   PetscFunctionReturn(PETSC_SUCCESS);
323 }
324 
325 static PetscErrorCode PCCompositeAddPC_Composite(PC pc, PC subpc)
326 {
327   PC_Composite    *jac;
328   PC_CompositeLink next, ilink;
329   PetscInt         cnt = 0;
330   const char      *prefix;
331   char             newprefix[20];
332 
333   PetscFunctionBegin;
334   PetscCall(PetscNew(&ilink));
335   ilink->next = NULL;
336   ilink->pc   = subpc;
337 
338   jac  = (PC_Composite *)pc->data;
339   next = jac->head;
340   if (!next) {
341     jac->head       = ilink;
342     ilink->previous = NULL;
343   } else {
344     cnt++;
345     while (next->next) {
346       next = next->next;
347       cnt++;
348     }
349     next->next      = ilink;
350     ilink->previous = next;
351   }
352   PetscCall(PCGetOptionsPrefix(pc, &prefix));
353   PetscCall(PCSetOptionsPrefix(subpc, prefix));
354   PetscCall(PetscSNPrintf(newprefix, 20, "sub_%d_", (int)cnt));
355   PetscCall(PCAppendOptionsPrefix(subpc, newprefix));
356   PetscCall(PetscObjectReference((PetscObject)subpc));
357   PetscFunctionReturn(PETSC_SUCCESS);
358 }
359 
360 static PetscErrorCode PCCompositeAddPCType_Composite(PC pc, PCType type)
361 {
362   PC subpc;
363 
364   PetscFunctionBegin;
365   PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), &subpc));
366   PetscCall(PetscObjectIncrementTabLevel((PetscObject)subpc, (PetscObject)pc, 1));
367   PetscCall(PCCompositeAddPC_Composite(pc, subpc));
368   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
369   PetscCall(PCSetType(subpc, type));
370   PetscCall(PCDestroy(&subpc));
371   PetscFunctionReturn(PETSC_SUCCESS);
372 }
373 
374 static PetscErrorCode PCCompositeGetNumberPC_Composite(PC pc, PetscInt *n)
375 {
376   PC_Composite    *jac;
377   PC_CompositeLink next;
378 
379   PetscFunctionBegin;
380   jac  = (PC_Composite *)pc->data;
381   next = jac->head;
382   *n   = 0;
383   while (next) {
384     next = next->next;
385     (*n)++;
386   }
387   PetscFunctionReturn(PETSC_SUCCESS);
388 }
389 
390 static PetscErrorCode PCCompositeGetPC_Composite(PC pc, PetscInt n, PC *subpc)
391 {
392   PC_Composite    *jac;
393   PC_CompositeLink next;
394   PetscInt         i;
395 
396   PetscFunctionBegin;
397   jac  = (PC_Composite *)pc->data;
398   next = jac->head;
399   for (i = 0; i < n; i++) {
400     PetscCheck(next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Not enough PCs in composite preconditioner");
401     next = next->next;
402   }
403   *subpc = next->pc;
404   PetscFunctionReturn(PETSC_SUCCESS);
405 }
406 
407 /*@
408   PCCompositeSetType - Sets the type of composite preconditioner.
409 
410   Logically Collective
411 
412   Input Parameters:
413 + pc   - the preconditioner context
414 - type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`
415 
416   Options Database Key:
417 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
418 
419   Level: advanced
420 
421 .seealso: `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`,
422           `PCCompositeGetType()`
423 @*/
424 PetscErrorCode PCCompositeSetType(PC pc, PCCompositeType type)
425 {
426   PetscFunctionBegin;
427   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
428   PetscValidLogicalCollectiveEnum(pc, type, 2);
429   PetscTryMethod(pc, "PCCompositeSetType_C", (PC, PCCompositeType), (pc, type));
430   PetscFunctionReturn(PETSC_SUCCESS);
431 }
432 
433 /*@
434   PCCompositeGetType - Gets the type of composite preconditioner.
435 
436   Logically Collective
437 
438   Input Parameter:
439 . pc - the preconditioner context
440 
441   Output Parameter:
442 . type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`
443 
444   Level: advanced
445 
446 .seealso: `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`,
447           `PCCompositeSetType()`
448 @*/
449 PetscErrorCode PCCompositeGetType(PC pc, PCCompositeType *type)
450 {
451   PetscFunctionBegin;
452   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
453   PetscUseMethod(pc, "PCCompositeGetType_C", (PC, PCCompositeType *), (pc, type));
454   PetscFunctionReturn(PETSC_SUCCESS);
455 }
456 
457 /*@
458   PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner, `PC_COMPOSITE_SPECIAL`,
459   for alphaI + R + S
460 
461   Logically Collective
462 
463   Input Parameters:
464 + pc    - the preconditioner context
465 - alpha - scale on identity
466 
467   Level: developer
468 
469 .seealso: `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`,
470           `PCCompositeSetType()`, `PCCompositeGetType()`
471 @*/
472 PetscErrorCode PCCompositeSpecialSetAlpha(PC pc, PetscScalar alpha)
473 {
474   PetscFunctionBegin;
475   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
476   PetscValidLogicalCollectiveScalar(pc, alpha, 2);
477   PetscTryMethod(pc, "PCCompositeSpecialSetAlpha_C", (PC, PetscScalar), (pc, alpha));
478   PetscFunctionReturn(PETSC_SUCCESS);
479 }
480 
481 /*@C
482   PCCompositeAddPCType - Adds another `PC` of the given type to the composite `PC`.
483 
484   Collective
485 
486   Input Parameters:
487 + pc   - the preconditioner context
488 - type - the type of the new preconditioner
489 
490   Level: intermediate
491 
492 .seealso: `PCCOMPOSITE`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()`
493 @*/
494 PetscErrorCode PCCompositeAddPCType(PC pc, PCType type)
495 {
496   PetscFunctionBegin;
497   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
498   PetscTryMethod(pc, "PCCompositeAddPCType_C", (PC, PCType), (pc, type));
499   PetscFunctionReturn(PETSC_SUCCESS);
500 }
501 
502 /*@
503   PCCompositeAddPC - Adds another `PC` to the composite `PC`.
504 
505   Collective
506 
507   Input Parameters:
508 + pc    - the preconditioner context
509 - subpc - the new preconditioner
510 
511   Level: intermediate
512 
513 .seealso: `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`
514 @*/
515 PetscErrorCode PCCompositeAddPC(PC pc, PC subpc)
516 {
517   PetscFunctionBegin;
518   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
519   PetscValidHeaderSpecific(subpc, PC_CLASSID, 2);
520   PetscTryMethod(pc, "PCCompositeAddPC_C", (PC, PC), (pc, subpc));
521   PetscFunctionReturn(PETSC_SUCCESS);
522 }
523 
524 /*@
525   PCCompositeGetNumberPC - Gets the number of `PC` objects in the composite `PC`.
526 
527   Not Collective
528 
529   Input Parameter:
530 . pc - the preconditioner context
531 
532   Output Parameter:
533 . num - the number of sub pcs
534 
535   Level: developer
536 
537 .seealso: `PCCOMPOSITE`, `PCCompositeGetPC()`, `PCCompositeAddPC()`, `PCCompositeAddPCType()`
538 @*/
539 PetscErrorCode PCCompositeGetNumberPC(PC pc, PetscInt *num)
540 {
541   PetscFunctionBegin;
542   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
543   PetscValidIntPointer(num, 2);
544   PetscUseMethod(pc, "PCCompositeGetNumberPC_C", (PC, PetscInt *), (pc, num));
545   PetscFunctionReturn(PETSC_SUCCESS);
546 }
547 
548 /*@
549   PCCompositeGetPC - Gets one of the `PC` objects in the composite `PC`.
550 
551   Not Collective
552 
553   Input Parameters:
554 + pc - the preconditioner context
555 - n  - the number of the pc requested
556 
557   Output Parameter:
558 . subpc - the PC requested
559 
560   Level: intermediate
561 
562   Note:
563   To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then
564   call `PCSetOperators()` on that `PC`.
565 
566 .seealso: `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`, `PCSetOperators()`
567 @*/
568 PetscErrorCode PCCompositeGetPC(PC pc, PetscInt n, PC *subpc)
569 {
570   PetscFunctionBegin;
571   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
572   PetscValidPointer(subpc, 3);
573   PetscUseMethod(pc, "PCCompositeGetPC_C", (PC, PetscInt, PC *), (pc, n, subpc));
574   PetscFunctionReturn(PETSC_SUCCESS);
575 }
576 
577 /*MC
578      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
579 
580    Options Database Keys:
581 +  -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
582 .  -pc_use_amat - activates `PCSetUseAmat()`
583 -  -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose
584 
585    Level: intermediate
586 
587    Notes:
588    To use a Krylov method inside the composite preconditioner, set the `PCType` of one or more
589    inner `PC`s to be `PCKSP`. Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
590    the incorrect answer) unless you use `KSPFGMRES` as the outer Krylov method
591 
592    To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then
593    call `PCSetOperators()` on that `PC`.
594 
595 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
596           `PCSHELL`, `PCKSP`, `PCCompositeSetType()`, `PCCompositeSpecialSetAlpha()`, `PCCompositeAddPCType()`,
597           `PCCompositeGetPC()`, `PCSetUseAmat()`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()`
598 M*/
599 
600 PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc)
601 {
602   PC_Composite *jac;
603 
604   PetscFunctionBegin;
605   PetscCall(PetscNew(&jac));
606 
607   pc->ops->apply           = PCApply_Composite_Additive;
608   pc->ops->applytranspose  = PCApplyTranspose_Composite_Additive;
609   pc->ops->setup           = PCSetUp_Composite;
610   pc->ops->setuponblocks   = PCSetUpOnBlocks_Composite;
611   pc->ops->reset           = PCReset_Composite;
612   pc->ops->destroy         = PCDestroy_Composite;
613   pc->ops->setfromoptions  = PCSetFromOptions_Composite;
614   pc->ops->view            = PCView_Composite;
615   pc->ops->applyrichardson = NULL;
616 
617   pc->data   = (void *)jac;
618   jac->type  = PC_COMPOSITE_ADDITIVE;
619   jac->work1 = NULL;
620   jac->work2 = NULL;
621   jac->head  = NULL;
622 
623   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", PCCompositeSetType_Composite));
624   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", PCCompositeGetType_Composite));
625   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", PCCompositeAddPCType_Composite));
626   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", PCCompositeAddPC_Composite));
627   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", PCCompositeGetNumberPC_Composite));
628   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", PCCompositeGetPC_Composite));
629   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", PCCompositeSpecialSetAlpha_Composite));
630   PetscFunctionReturn(PETSC_SUCCESS);
631 }
632