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