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