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