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