Actual source code: composite.c
1: #define PETSCKSP_DLL
3: /*
4: Defines a preconditioner that can consist of a collection of PCs
5: */
6: #include private/pcimpl.h
7: #include petscksp.h
9: typedef struct _PC_CompositeLink *PC_CompositeLink;
10: struct _PC_CompositeLink {
11: PC pc;
12: PC_CompositeLink next;
13: PC_CompositeLink previous;
14: };
15:
16: typedef struct {
17: PC_CompositeLink head;
18: PCCompositeType type;
19: Vec work1;
20: Vec work2;
21: PetscScalar alpha;
22: PetscTruth use_true_matrix;
23: } PC_Composite;
27: static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
28: {
29: PetscErrorCode ierr;
30: PC_Composite *jac = (PC_Composite*)pc->data;
31: PC_CompositeLink next = jac->head;
32: Mat mat = pc->pmat;
35: if (!next) {
36: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
37: }
38: if (next->next && !jac->work2) { /* allocate second work vector */
39: VecDuplicate(jac->work1,&jac->work2);
40: }
41: PCApply(next->pc,x,y);
42: if (jac->use_true_matrix) mat = pc->mat;
43: while (next->next) {
44: next = next->next;
45: MatMult(mat,y,jac->work1);
46: VecWAXPY(jac->work2,-1.0,jac->work1,x);
47: VecSet(jac->work1,0.0); /* zero since some PC's may not set all entries in the result */
48: PCApply(next->pc,jac->work2,jac->work1);
49: VecAXPY(y,1.0,jac->work1);
50: }
51: if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
52: while (next->previous) {
53: next = next->previous;
54: MatMult(mat,y,jac->work1);
55: VecWAXPY(jac->work2,-1.0,jac->work1,x);
56: VecSet(jac->work1,0.0); /* zero since some PC's may not set all entries in the result */
57: PCApply(next->pc,jac->work2,jac->work1);
58: VecAXPY(y,1.0,jac->work1);
59: }
60: }
61: return(0);
62: }
64: /*
65: This is very special for a matrix of the form alpha I + R + S
66: where first preconditioner is built from alpha I + S and second from
67: alpha I + R
68: */
71: static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
72: {
73: PetscErrorCode ierr;
74: PC_Composite *jac = (PC_Composite*)pc->data;
75: PC_CompositeLink next = jac->head;
78: if (!next) {
79: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
80: }
81: if (!next->next || next->next->next) {
82: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
83: }
85: PCApply(next->pc,x,jac->work1);
86: PCApply(next->next->pc,jac->work1,y);
87: return(0);
88: }
92: static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
93: {
94: PetscErrorCode ierr;
95: PC_Composite *jac = (PC_Composite*)pc->data;
96: PC_CompositeLink next = jac->head;
99: if (!next) {
100: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
101: }
102: PCApply(next->pc,x,y);
103: while (next->next) {
104: next = next->next;
105: VecSet(jac->work1,0.0); /* zero since some PC's may not set all entries in the result */
106: PCApply(next->pc,x,jac->work1);
107: VecAXPY(y,1.0,jac->work1);
108: }
109: return(0);
110: }
114: static PetscErrorCode PCSetUp_Composite(PC pc)
115: {
116: PetscErrorCode ierr;
117: PC_Composite *jac = (PC_Composite*)pc->data;
118: PC_CompositeLink next = jac->head;
121: if (!jac->work1) {
122: MatGetVecs(pc->pmat,&jac->work1,0);
123: }
124: while (next) {
125: PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);
126: next = next->next;
127: }
128: return(0);
129: }
133: static PetscErrorCode PCDestroy_Composite(PC pc)
134: {
135: PC_Composite *jac = (PC_Composite*)pc->data;
136: PetscErrorCode ierr;
137: PC_CompositeLink next = jac->head;
140: while (next) {
141: PCDestroy(next->pc);
142: next = next->next;
143: }
145: if (jac->work1) {VecDestroy(jac->work1);}
146: if (jac->work2) {VecDestroy(jac->work2);}
147: PetscFree(jac);
148: return(0);
149: }
153: static PetscErrorCode PCSetFromOptions_Composite(PC pc)
154: {
155: PC_Composite *jac = (PC_Composite*)pc->data;
156: PetscErrorCode ierr;
157: PetscInt nmax = 8,i;
158: PC_CompositeLink next;
159: char *pcs[8];
160: PetscTruth flg;
163: PetscOptionsHead("Composite preconditioner options");
164: PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);
165: if (flg) {
166: PCCompositeSetType(pc,jac->type);
167: }
168: flg = PETSC_FALSE;
169: PetscOptionsTruth("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",flg,&flg,PETSC_NULL);
170: if (flg) {
171: PCCompositeSetUseTrue(pc);
172: }
173: PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);
174: if (flg) {
175: for (i=0; i<nmax; i++) {
176: PCCompositeAddPC(pc,pcs[i]);
177: }
178: }
179: PetscOptionsTail();
181: next = jac->head;
182: while (next) {
183: PCSetFromOptions(next->pc);
184: next = next->next;
185: }
186: return(0);
187: }
191: static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
192: {
193: PC_Composite *jac = (PC_Composite*)pc->data;
194: PetscErrorCode ierr;
195: PC_CompositeLink next = jac->head;
196: PetscTruth iascii;
199: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
200: if (iascii) {
201: PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);
202: PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");
203: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
204: } else {
205: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
206: }
207: if (iascii) {
208: PetscViewerASCIIPushTab(viewer);
209: }
210: while (next) {
211: PCView(next->pc,viewer);
212: next = next->next;
213: }
214: if (iascii) {
215: PetscViewerASCIIPopTab(viewer);
216: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
217: }
218: return(0);
219: }
221: /* ------------------------------------------------------------------------------*/
226: PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
227: {
228: PC_Composite *jac = (PC_Composite*)pc->data;
230: jac->alpha = alpha;
231: return(0);
232: }
238: PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type)
239: {
241: if (type == PC_COMPOSITE_ADDITIVE) {
242: pc->ops->apply = PCApply_Composite_Additive;
243: } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
244: pc->ops->apply = PCApply_Composite_Multiplicative;
245: } else if (type == PC_COMPOSITE_SPECIAL) {
246: pc->ops->apply = PCApply_Composite_Special;
247: } else {
248: SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
249: }
250: return(0);
251: }
257: PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type)
258: {
259: PC_Composite *jac;
260: PC_CompositeLink next,ilink;
261: PetscErrorCode ierr;
262: PetscInt cnt = 0;
263: const char *prefix;
264: char newprefix[8];
267: PetscNewLog(pc,struct _PC_CompositeLink,&ilink);
268: ilink->next = 0;
269: PCCreate(((PetscObject)pc)->comm,&ilink->pc);
271: jac = (PC_Composite*)pc->data;
272: next = jac->head;
273: if (!next) {
274: jac->head = ilink;
275: ilink->previous = PETSC_NULL;
276: } else {
277: cnt++;
278: while (next->next) {
279: next = next->next;
280: cnt++;
281: }
282: next->next = ilink;
283: ilink->previous = next;
284: }
285: PCGetOptionsPrefix(pc,&prefix);
286: PCSetOptionsPrefix(ilink->pc,prefix);
287: sprintf(newprefix,"sub_%d_",(int)cnt);
288: PCAppendOptionsPrefix(ilink->pc,newprefix);
289: /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
290: PCSetType(ilink->pc,type);
291: return(0);
292: }
298: PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
299: {
300: PC_Composite *jac;
301: PC_CompositeLink next;
302: PetscInt i;
305: jac = (PC_Composite*)pc->data;
306: next = jac->head;
307: for (i=0; i<n; i++) {
308: if (!next->next) {
309: SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
310: }
311: next = next->next;
312: }
313: *subpc = next->pc;
314: return(0);
315: }
321: PetscErrorCode PCCompositeSetUseTrue_Composite(PC pc)
322: {
323: PC_Composite *jac;
326: jac = (PC_Composite*)pc->data;
327: jac->use_true_matrix = PETSC_TRUE;
328: return(0);
329: }
332: /* -------------------------------------------------------------------------------- */
335: /*@
336: PCCompositeSetType - Sets the type of composite preconditioner.
337:
338: Collective on PC
340: Input Parameter:
341: + pc - the preconditioner context
342: - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
344: Options Database Key:
345: . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
347: Level: Developer
349: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
350: @*/
351: PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type)
352: {
353: PetscErrorCode ierr,(*f)(PC,PCCompositeType);
357: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);
358: if (f) {
359: (*f)(pc,type);
360: }
361: return(0);
362: }
366: /*@
367: PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
368: for alphaI + R + S
369:
370: Collective on PC
372: Input Parameter:
373: + pc - the preconditioner context
374: - alpha - scale on identity
376: Level: Developer
378: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
379: @*/
380: PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
381: {
382: PetscErrorCode ierr,(*f)(PC,PetscScalar);
386: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);
387: if (f) {
388: (*f)(pc,alpha);
389: }
390: return(0);
391: }
395: /*@C
396: PCCompositeAddPC - Adds another PC to the composite PC.
397:
398: Collective on PC
400: Input Parameters:
401: + pc - the preconditioner context
402: - type - the type of the new preconditioner
404: Level: Developer
406: .keywords: PC, composite preconditioner, add
407: @*/
408: PetscErrorCode PCCompositeAddPC(PC pc,PCType type)
409: {
410: PetscErrorCode ierr,(*f)(PC,PCType);
414: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);
415: if (f) {
416: (*f)(pc,type);
417: }
418: return(0);
419: }
423: /*@
424: PCCompositeGetPC - Gets one of the PC objects in the composite PC.
425:
426: Not Collective
428: Input Parameter:
429: + pc - the preconditioner context
430: - n - the number of the pc requested
432: Output Parameters:
433: . subpc - the PC requested
435: Level: Developer
437: .keywords: PC, get, composite preconditioner, sub preconditioner
439: .seealso: PCCompositeAddPC()
440: @*/
441: PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
442: {
443: PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
448: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);
449: if (f) {
450: (*f)(pc,n,subpc);
451: } else {
452: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
453: }
454: return(0);
455: }
459: /*@
460: PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
461: the matrix used to define the preconditioner) is used to compute
462: the residual when the multiplicative scheme is used.
464: Collective on PC
466: Input Parameters:
467: . pc - the preconditioner context
469: Options Database Key:
470: . -pc_composite_true - Activates PCCompositeSetUseTrue()
472: Note:
473: For the common case in which the preconditioning and linear
474: system matrices are identical, this routine is unnecessary.
476: Level: Developer
478: .keywords: PC, composite preconditioner, set, true, flag
480: .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
481: @*/
482: PetscErrorCode PCCompositeSetUseTrue(PC pc)
483: {
484: PetscErrorCode ierr,(*f)(PC);
488: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);
489: if (f) {
490: (*f)(pc);
491: }
492: return(0);
493: }
495: /* -------------------------------------------------------------------------------------------*/
497: /*MC
498: PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
500: Options Database Keys:
501: + -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
502: . -pc_composite_true - Activates PCCompositeSetUseTrue()
503: - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose
505: Level: intermediate
507: Concepts: composing solvers
509: Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
510: inner PCs to be PCKSP.
511: Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
512: the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
515: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
516: PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
517: PCCompositeGetPC(), PCCompositeSetUseTrue()
519: M*/
524: PetscErrorCode PCCreate_Composite(PC pc)
525: {
527: PC_Composite *jac;
530: PetscNewLog(pc,PC_Composite,&jac);
531: pc->ops->apply = PCApply_Composite_Additive;
532: pc->ops->setup = PCSetUp_Composite;
533: pc->ops->destroy = PCDestroy_Composite;
534: pc->ops->setfromoptions = PCSetFromOptions_Composite;
535: pc->ops->view = PCView_Composite;
536: pc->ops->applyrichardson = 0;
538: pc->data = (void*)jac;
539: jac->type = PC_COMPOSITE_ADDITIVE;
540: jac->work1 = 0;
541: jac->work2 = 0;
542: jac->head = 0;
544: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
545: PCCompositeSetType_Composite);
546: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
547: PCCompositeAddPC_Composite);
548: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
549: PCCompositeGetPC_Composite);
550: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
551: PCCompositeSetUseTrue_Composite);
552: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
553: PCCompositeSpecialSetAlpha_Composite);
555: return(0);
556: }