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