Actual source code: composite.c
1: /*$Id: composite.c,v 1.45 2001/08/07 03:03:39 balay Exp $*/
2: /*
3: Defines a preconditioner that can consist of a collection of PCs
4: */
5: #include src/ksp/pc/pcimpl.h
6: #include petscksp.h
8: typedef struct _PC_CompositeLink *PC_CompositeLink;
9: struct _PC_CompositeLink {
10: PC pc;
11: PC_CompositeLink next;
12: };
13:
14: typedef struct {
15: PC_CompositeLink head;
16: PCCompositeType type;
17: Vec work1;
18: Vec work2;
19: PetscScalar alpha;
20: PetscTruth use_true_matrix;
21: } PC_Composite;
25: static int PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
26: {
27: int ierr;
28: PC_Composite *jac = (PC_Composite*)pc->data;
29: PC_CompositeLink next = jac->head;
30: PetscScalar one = 1.0,mone = -1.0;
31: Mat mat = pc->pmat;
34: if (!next) {
35: SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()");
36: }
37: if (next->next && !jac->work2) { /* allocate second work vector */
38: VecDuplicate(jac->work1,&jac->work2);
39: }
40: PCApply(next->pc,x,y,PC_LEFT);
41: if (jac->use_true_matrix) mat = pc->mat;
42: while (next->next) {
43: next = next->next;
44: MatMult(mat,y,jac->work1);
45: VecWAXPY(&mone,jac->work1,x,jac->work2);
46: PCApply(next->pc,jac->work2,jac->work1,PC_LEFT);
47: VecAXPY(&one,jac->work1,y);
48: }
50: return(0);
51: }
53: /*
54: This is very special for a matrix of the form alpha I + R + S
55: where first preconditioner is built from alpha I + S and second from
56: alpha I + R
57: */
60: static int PCApply_Composite_Special(PC pc,Vec x,Vec y)
61: {
62: int ierr;
63: PC_Composite *jac = (PC_Composite*)pc->data;
64: PC_CompositeLink next = jac->head;
67: if (!next) {
68: SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()");
69: }
70: if (!next->next || next->next->next) {
71: SETERRQ(1,"Special composite preconditioners requires exactly two PCs");
72: }
74: PCApply(next->pc,x,jac->work1,PC_LEFT);
75: PCApply(next->next->pc,jac->work1,y,PC_LEFT);
76: return(0);
77: }
81: static int PCApply_Composite_Additive(PC pc,Vec x,Vec y)
82: {
83: int ierr;
84: PC_Composite *jac = (PC_Composite*)pc->data;
85: PC_CompositeLink next = jac->head;
86: PetscScalar one = 1.0;
89: if (!next) {
90: SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()");
91: }
92: PCApply(next->pc,x,y,PC_LEFT);
93: while (next->next) {
94: next = next->next;
95: PCApply(next->pc,x,jac->work1,PC_LEFT);
96: VecAXPY(&one,jac->work1,y);
97: }
98: return(0);
99: }
103: static int PCSetUp_Composite(PC pc)
104: {
105: int ierr;
106: PC_Composite *jac = (PC_Composite*)pc->data;
107: PC_CompositeLink next = jac->head;
110: if (!jac->work1) {
111: VecDuplicate(pc->vec,&jac->work1);
112: }
113: while (next) {
114: PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);
115: PCSetVector(next->pc,jac->work1);
116: next = next->next;
117: }
119: return(0);
120: }
124: static int PCDestroy_Composite(PC pc)
125: {
126: PC_Composite *jac = (PC_Composite*)pc->data;
127: int ierr;
128: PC_CompositeLink next = jac->head;
131: while (next) {
132: PCDestroy(next->pc);
133: next = next->next;
134: }
136: if (jac->work1) {VecDestroy(jac->work1);}
137: if (jac->work2) {VecDestroy(jac->work2);}
138: PetscFree(jac);
139: return(0);
140: }
144: static int PCSetFromOptions_Composite(PC pc)
145: {
146: PC_Composite *jac = (PC_Composite*)pc->data;
147: int ierr,nmax = 8,i,indx;
148: PC_CompositeLink next;
149: char *pcs[8];
150: const char *types[] = {"multiplicative","additive","special"};
151: PetscTruth flg;
154: PetscOptionsHead("Composite preconditioner options");
155: PetscOptionsEList("-pc_composite_type","Type of composition","PCCompositeSetType",types,3,types[0],&indx,&flg);
156: if (flg) {
157: PCCompositeSetType(pc,(PCCompositeType)indx);
158: }
159: PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);
160: if (flg) {
161: PCCompositeSetUseTrue(pc);
162: }
163: PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);
164: if (flg) {
165: for (i=0; i<nmax; i++) {
166: PCCompositeAddPC(pc,pcs[i]);
167: }
168: }
169: PetscOptionsTail();
171: next = jac->head;
172: while (next) {
173: PCSetFromOptions(next->pc);
174: next = next->next;
175: }
176: return(0);
177: }
181: static int PCView_Composite(PC pc,PetscViewer viewer)
182: {
183: PC_Composite *jac = (PC_Composite*)pc->data;
184: int ierr;
185: PC_CompositeLink next = jac->head;
186: PetscTruth isascii;
189: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
190: if (isascii) {
191: PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");
192: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
193: } else {
194: SETERRQ1(1,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
195: }
196: if (isascii) {
197: PetscViewerASCIIPushTab(viewer);
198: }
199: while (next) {
200: PCView(next->pc,viewer);
201: next = next->next;
202: }
203: if (isascii) {
204: PetscViewerASCIIPopTab(viewer);
205: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
206: }
207: return(0);
208: }
210: /* ------------------------------------------------------------------------------*/
212: EXTERN_C_BEGIN
215: int PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
216: {
217: PC_Composite *jac = (PC_Composite*)pc->data;
219: jac->alpha = alpha;
220: return(0);
221: }
222: EXTERN_C_END
224: EXTERN_C_BEGIN
227: int PCCompositeSetType_Composite(PC pc,PCCompositeType type)
228: {
230: if (type == PC_COMPOSITE_ADDITIVE) {
231: pc->ops->apply = PCApply_Composite_Additive;
232: } else if (type == PC_COMPOSITE_MULTIPLICATIVE) {
233: pc->ops->apply = PCApply_Composite_Multiplicative;
234: } else if (type == PC_COMPOSITE_SPECIAL) {
235: pc->ops->apply = PCApply_Composite_Special;
236: } else {
237: SETERRQ(1,"Unkown composite preconditioner type");
238: }
239: return(0);
240: }
241: EXTERN_C_END
243: EXTERN_C_BEGIN
246: int PCCompositeAddPC_Composite(PC pc,PCType type)
247: {
248: PC_Composite *jac;
249: PC_CompositeLink next,link;
250: int ierr,cnt = 0;
251: char *prefix,newprefix[8];
254: PetscNew(struct _PC_CompositeLink,&link);
255: link->next = 0;
256: PCCreate(pc->comm,&link->pc);
258: jac = (PC_Composite*)pc->data;
259: next = jac->head;
260: if (!next) {
261: jac->head = link;
262: } else {
263: cnt++;
264: while (next->next) {
265: next = next->next;
266: cnt++;
267: }
268: next->next = link;
269: }
270: PCGetOptionsPrefix(pc,&prefix);
271: PCSetOptionsPrefix(link->pc,prefix);
272: sprintf(newprefix,"sub_%d_",cnt);
273: PCAppendOptionsPrefix(link->pc,newprefix);
274: /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
275: PCSetType(link->pc,type);
277: return(0);
278: }
279: EXTERN_C_END
281: EXTERN_C_BEGIN
284: int PCCompositeGetPC_Composite(PC pc,int n,PC *subpc)
285: {
286: PC_Composite *jac;
287: PC_CompositeLink next;
288: int i;
291: jac = (PC_Composite*)pc->data;
292: next = jac->head;
293: for (i=0; i<n; i++) {
294: if (!next->next) {
295: SETERRQ(1,"Not enough PCs in composite preconditioner");
296: }
297: next = next->next;
298: }
299: *subpc = next->pc;
300: return(0);
301: }
302: EXTERN_C_END
304: EXTERN_C_BEGIN
307: int PCCompositeSetUseTrue_Composite(PC pc)
308: {
309: PC_Composite *jac;
312: jac = (PC_Composite*)pc->data;
313: jac->use_true_matrix = PETSC_TRUE;
314: return(0);
315: }
316: EXTERN_C_END
318: /* -------------------------------------------------------------------------------- */
321: /*@C
322: PCCompositeSetType - Sets the type of composite preconditioner.
323:
324: Collective on PC
326: Input Parameter:
327: . pc - the preconditioner context
328: . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
330: Options Database Key:
331: . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
333: Level: Developer
335: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
336: @*/
337: int PCCompositeSetType(PC pc,PCCompositeType type)
338: {
339: int ierr,(*f)(PC,PCCompositeType);
343: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);
344: if (f) {
345: (*f)(pc,type);
346: }
347: return(0);
348: }
352: /*@C
353: PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
354: for alphaI + R + S
355:
356: Collective on PC
358: Input Parameter:
359: + pc - the preconditioner context
360: - alpha - scale on identity
362: Level: Developer
364: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
365: @*/
366: int PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
367: {
368: int ierr,(*f)(PC,PetscScalar);
372: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);
373: if (f) {
374: (*f)(pc,alpha);
375: }
376: return(0);
377: }
381: /*@C
382: PCCompositeAddPC - Adds another PC to the composite PC.
383:
384: Collective on PC
386: Input Parameters:
387: . pc - the preconditioner context
388: . type - the type of the new preconditioner
390: Level: Developer
392: .keywords: PC, composite preconditioner, add
393: @*/
394: int PCCompositeAddPC(PC pc,PCType type)
395: {
396: int ierr,(*f)(PC,PCType);
400: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);
401: if (f) {
402: (*f)(pc,type);
403: }
404: return(0);
405: }
409: /*@C
410: PCCompositeGetPC - Gets one of the PC objects in the composite PC.
411:
412: Not Collective
414: Input Parameter:
415: . pc - the preconditioner context
416: . n - the number of the pc requested
418: Output Parameters:
419: . subpc - the PC requested
421: Level: Developer
423: .keywords: PC, get, composite preconditioner, sub preconditioner
425: .seealso: PCCompositeAddPC()
426: @*/
427: int PCCompositeGetPC(PC pc,int n,PC *subpc)
428: {
429: int ierr,(*f)(PC,int,PC *);
434: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);
435: if (f) {
436: (*f)(pc,n,subpc);
437: } else {
438: SETERRQ(1,"Cannot get pc, not composite type");
439: }
440: return(0);
441: }
445: /*@
446: PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
447: the matrix used to define the preconditioner) is used to compute
448: the residual when the multiplicative scheme is used.
450: Collective on PC
452: Input Parameters:
453: . pc - the preconditioner context
455: Options Database Key:
456: . -pc_composite_true - Activates PCCompositeSetUseTrue()
458: Note:
459: For the common case in which the preconditioning and linear
460: system matrices are identical, this routine is unnecessary.
462: Level: Developer
464: .keywords: PC, composite preconditioner, set, true, flag
466: .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
467: @*/
468: int PCCompositeSetUseTrue(PC pc)
469: {
470: int ierr,(*f)(PC);
474: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);
475: if (f) {
476: (*f)(pc);
477: }
478: return(0);
479: }
481: /* -------------------------------------------------------------------------------------------*/
483: /*MC
484: PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
486: Options Database Keys:
487: . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
488: . -pc_composite_true - Activates PCCompositeSetUseTrue()
490: Level: intermediate
492: Concepts: composing solvers
494: Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
495: inner PCs to be PCKSP.
496: Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
497: the incorrect answer) unless you use KSPFGMRES as the other Krylov method
500: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
501: PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
502: PCCompositeGetPC(), PCCompositeSetUseTrue()
504: M*/
506: EXTERN_C_BEGIN
509: int PCCreate_Composite(PC pc)
510: {
511: int ierr;
512: PC_Composite *jac;
515: PetscNew(PC_Composite,&jac);
516: PetscLogObjectMemory(pc,sizeof(PC_Composite));
517: pc->ops->apply = PCApply_Composite_Additive;
518: pc->ops->setup = PCSetUp_Composite;
519: pc->ops->destroy = PCDestroy_Composite;
520: pc->ops->setfromoptions = PCSetFromOptions_Composite;
521: pc->ops->view = PCView_Composite;
522: pc->ops->applyrichardson = 0;
524: pc->data = (void*)jac;
525: jac->type = PC_COMPOSITE_ADDITIVE;
526: jac->work1 = 0;
527: jac->work2 = 0;
528: jac->head = 0;
530: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
531: PCCompositeSetType_Composite);
532: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
533: PCCompositeAddPC_Composite);
534: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
535: PCCompositeGetPC_Composite);
536: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
537: PCCompositeSetUseTrue_Composite);
538: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
539: PCCompositeSpecialSetAlpha_Composite);
541: return(0);
542: }
543: EXTERN_C_END