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