Actual source code: composite.c
1: /*
2: Defines a preconditioner that can consist of a collection of PCs
3: */
4: #include src/ksp/pc/pcimpl.h
5: #include petscksp.h
7: typedef struct _PC_CompositeLink *PC_CompositeLink;
8: struct _PC_CompositeLink {
9: PC pc;
10: PC_CompositeLink next;
11: };
12:
13: typedef struct {
14: PC_CompositeLink head;
15: PCCompositeType type;
16: Vec work1;
17: Vec work2;
18: PetscScalar alpha;
19: PetscTruth use_true_matrix;
20: } PC_Composite;
24: static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
25: {
26: PetscErrorCode ierr;
27: PC_Composite *jac = (PC_Composite*)pc->data;
28: PC_CompositeLink next = jac->head;
29: PetscScalar one = 1.0,mone = -1.0;
30: Mat mat = pc->pmat;
33: if (!next) {
34: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
35: }
36: if (next->next && !jac->work2) { /* allocate second work vector */
37: VecDuplicate(jac->work1,&jac->work2);
38: }
39: PCApply(next->pc,x,y);
40: if (jac->use_true_matrix) mat = pc->mat;
41: while (next->next) {
42: next = next->next;
43: MatMult(mat,y,jac->work1);
44: VecWAXPY(&mone,jac->work1,x,jac->work2);
45: PCApply(next->pc,jac->work2,jac->work1);
46: VecAXPY(&one,jac->work1,y);
47: }
49: return(0);
50: }
52: /*
53: This is very special for a matrix of the form alpha I + R + S
54: where first preconditioner is built from alpha I + S and second from
55: alpha I + R
56: */
59: static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
60: {
61: PetscErrorCode ierr;
62: PC_Composite *jac = (PC_Composite*)pc->data;
63: PC_CompositeLink next = jac->head;
66: if (!next) {
67: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
68: }
69: if (!next->next || next->next->next) {
70: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
71: }
73: PCApply(next->pc,x,jac->work1);
74: PCApply(next->next->pc,jac->work1,y);
75: return(0);
76: }
80: static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
81: {
82: PetscErrorCode ierr;
83: PC_Composite *jac = (PC_Composite*)pc->data;
84: PC_CompositeLink next = jac->head;
85: PetscScalar one = 1.0;
88: if (!next) {
89: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
90: }
91: PCApply(next->pc,x,y);
92: while (next->next) {
93: next = next->next;
94: PCApply(next->pc,x,jac->work1);
95: VecAXPY(&one,jac->work1,y);
96: }
97: return(0);
98: }
102: static PetscErrorCode PCSetUp_Composite(PC pc)
103: {
104: PetscErrorCode ierr;
105: PC_Composite *jac = (PC_Composite*)pc->data;
106: PC_CompositeLink next = jac->head;
109: if (!jac->work1) {
110: MatGetVecs(pc->pmat,&jac->work1,0);
111: }
112: while (next) {
113: PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);
114: next = next->next;
115: }
117: return(0);
118: }
122: static PetscErrorCode PCDestroy_Composite(PC pc)
123: {
124: PC_Composite *jac = (PC_Composite*)pc->data;
125: PetscErrorCode ierr;
126: PC_CompositeLink next = jac->head;
129: while (next) {
130: PCDestroy(next->pc);
131: next = next->next;
132: }
134: if (jac->work1) {VecDestroy(jac->work1);}
135: if (jac->work2) {VecDestroy(jac->work2);}
136: PetscFree(jac);
137: return(0);
138: }
142: static PetscErrorCode PCSetFromOptions_Composite(PC pc)
143: {
144: PC_Composite *jac = (PC_Composite*)pc->data;
145: PetscErrorCode ierr;
146: PetscInt nmax = 8,i,indx;
147: PC_CompositeLink next;
148: char *pcs[8];
149: const char *types[] = {"additive","multiplicative","special"};
150: PetscTruth flg;
153: PetscOptionsHead("Composite preconditioner options");
154: PetscOptionsEList("-pc_composite_type","Type of composition","PCCompositeSetType",types,3,types[0],&indx,&flg);
155: if (flg) {
156: PCCompositeSetType(pc,(PCCompositeType)indx);
157: }
158: PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);
159: if (flg) {
160: PCCompositeSetUseTrue(pc);
161: }
162: PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);
163: if (flg) {
164: for (i=0; i<nmax; i++) {
165: PCCompositeAddPC(pc,pcs[i]);
166: }
167: }
168: PetscOptionsTail();
170: next = jac->head;
171: while (next) {
172: PCSetFromOptions(next->pc);
173: next = next->next;
174: }
175: return(0);
176: }
180: static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
181: {
182: PC_Composite *jac = (PC_Composite*)pc->data;
183: PetscErrorCode ierr;
184: PC_CompositeLink next = jac->head;
185: PetscTruth iascii;
188: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
189: if (iascii) {
190: PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");
191: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
192: } else {
193: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
194: }
195: if (iascii) {
196: PetscViewerASCIIPushTab(viewer);
197: }
198: while (next) {
199: PCView(next->pc,viewer);
200: next = next->next;
201: }
202: if (iascii) {
203: PetscViewerASCIIPopTab(viewer);
204: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
205: }
206: return(0);
207: }
209: /* ------------------------------------------------------------------------------*/
214: PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
215: {
216: PC_Composite *jac = (PC_Composite*)pc->data;
218: jac->alpha = alpha;
219: return(0);
220: }
226: PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type)
227: {
229: if (type == PC_COMPOSITE_ADDITIVE) {
230: pc->ops->apply = PCApply_Composite_Additive;
231: } else if (type == PC_COMPOSITE_MULTIPLICATIVE) {
232: pc->ops->apply = PCApply_Composite_Multiplicative;
233: } else if (type == PC_COMPOSITE_SPECIAL) {
234: pc->ops->apply = PCApply_Composite_Special;
235: } else {
236: SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
237: }
238: return(0);
239: }
245: PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type)
246: {
247: PC_Composite *jac;
248: PC_CompositeLink next,link;
249: PetscErrorCode ierr;
250: PetscInt 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_",(int)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: }
284: PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
285: {
286: PC_Composite *jac;
287: PC_CompositeLink next;
288: PetscInt i;
291: jac = (PC_Composite*)pc->data;
292: next = jac->head;
293: for (i=0; i<n; i++) {
294: if (!next->next) {
295: SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
296: }
297: next = next->next;
298: }
299: *subpc = next->pc;
300: return(0);
301: }
307: PetscErrorCode 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: }
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: PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type)
338: {
339: PetscErrorCode 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: PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
367: {
368: PetscErrorCode 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: PetscErrorCode PCCompositeAddPC(PC pc,PCType type)
395: {
396: PetscErrorCode 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: PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
428: {
429: PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
434: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);
435: if (f) {
436: (*f)(pc,n,subpc);
437: } else {
438: SETERRQ(PETSC_ERR_ARG_WRONG,"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: PetscErrorCode PCCompositeSetUseTrue(PC pc)
469: {
470: PetscErrorCode 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 outter Krylov method
500: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
501: PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
502: PCCompositeGetPC(), PCCompositeSetUseTrue()
504: M*/
509: PetscErrorCode PCCreate_Composite(PC pc)
510: {
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: }