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: };
14:
15: typedef struct {
16: PC_CompositeLink head;
17: PCCompositeType type;
18: Vec work1;
19: Vec work2;
20: PetscScalar alpha;
21: PetscTruth use_true_matrix;
22: } PC_Composite;
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;
34: if (!next) {
35: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"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);
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(jac->work2,-1.0,jac->work1,x);
46: PCApply(next->pc,jac->work2,jac->work1);
47: VecAXPY(y,1.0,jac->work1);
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 PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
61: {
62: PetscErrorCode ierr;
63: PC_Composite *jac = (PC_Composite*)pc->data;
64: PC_CompositeLink next = jac->head;
67: if (!next) {
68: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
69: }
70: if (!next->next || next->next->next) {
71: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
72: }
74: PCApply(next->pc,x,jac->work1);
75: PCApply(next->next->pc,jac->work1,y);
76: return(0);
77: }
81: static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
82: {
83: PetscErrorCode ierr;
84: PC_Composite *jac = (PC_Composite*)pc->data;
85: PC_CompositeLink next = jac->head;
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(y,1.0,jac->work1);
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;
147: PC_CompositeLink next;
148: char *pcs[8];
149: PetscTruth flg;
152: PetscOptionsHead("Composite preconditioner options");
153: PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);
154: PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);
155: if (flg) {
156: PCCompositeSetUseTrue(pc);
157: }
158: PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);
159: if (flg) {
160: for (i=0; i<nmax; i++) {
161: PCCompositeAddPC(pc,pcs[i]);
162: }
163: }
164: PetscOptionsTail();
166: next = jac->head;
167: while (next) {
168: PCSetFromOptions(next->pc);
169: next = next->next;
170: }
171: return(0);
172: }
176: static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
177: {
178: PC_Composite *jac = (PC_Composite*)pc->data;
179: PetscErrorCode ierr;
180: PC_CompositeLink next = jac->head;
181: PetscTruth iascii;
184: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
185: if (iascii) {
186: PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);
187: PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");
188: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
189: } else {
190: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
191: }
192: if (iascii) {
193: PetscViewerASCIIPushTab(viewer);
194: }
195: while (next) {
196: PCView(next->pc,viewer);
197: next = next->next;
198: }
199: if (iascii) {
200: PetscViewerASCIIPopTab(viewer);
201: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
202: }
203: return(0);
204: }
206: /* ------------------------------------------------------------------------------*/
211: PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
212: {
213: PC_Composite *jac = (PC_Composite*)pc->data;
215: jac->alpha = alpha;
216: return(0);
217: }
223: PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type)
224: {
226: if (type == PC_COMPOSITE_ADDITIVE) {
227: pc->ops->apply = PCApply_Composite_Additive;
228: } else if (type == PC_COMPOSITE_MULTIPLICATIVE) {
229: pc->ops->apply = PCApply_Composite_Multiplicative;
230: } else if (type == PC_COMPOSITE_SPECIAL) {
231: pc->ops->apply = PCApply_Composite_Special;
232: } else {
233: SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
234: }
235: return(0);
236: }
242: PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type)
243: {
244: PC_Composite *jac;
245: PC_CompositeLink next,ilink;
246: PetscErrorCode ierr;
247: PetscInt cnt = 0;
248: const char *prefix;
249: char newprefix[8];
252: PetscNew(struct _PC_CompositeLink,&ilink);
253: ilink->next = 0;
254: PCCreate(pc->comm,&ilink->pc);
256: jac = (PC_Composite*)pc->data;
257: next = jac->head;
258: if (!next) {
259: jac->head = ilink;
260: } else {
261: cnt++;
262: while (next->next) {
263: next = next->next;
264: cnt++;
265: }
266: next->next = ilink;
267: }
268: PCGetOptionsPrefix(pc,&prefix);
269: PCSetOptionsPrefix(ilink->pc,prefix);
270: sprintf(newprefix,"sub_%d_",(int)cnt);
271: PCAppendOptionsPrefix(ilink->pc,newprefix);
272: /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
273: PCSetType(ilink->pc,type);
275: return(0);
276: }
282: PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
283: {
284: PC_Composite *jac;
285: PC_CompositeLink next;
286: PetscInt i;
289: jac = (PC_Composite*)pc->data;
290: next = jac->head;
291: for (i=0; i<n; i++) {
292: if (!next->next) {
293: SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
294: }
295: next = next->next;
296: }
297: *subpc = next->pc;
298: return(0);
299: }
305: PetscErrorCode PCCompositeSetUseTrue_Composite(PC pc)
306: {
307: PC_Composite *jac;
310: jac = (PC_Composite*)pc->data;
311: jac->use_true_matrix = PETSC_TRUE;
312: return(0);
313: }
316: /* -------------------------------------------------------------------------------- */
319: /*@
320: PCCompositeSetType - Sets the type of composite preconditioner.
321:
322: Collective on PC
324: Input Parameter:
325: + pc - the preconditioner context
326: - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
328: Options Database Key:
329: . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
331: Level: Developer
333: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
334: @*/
335: PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type)
336: {
337: PetscErrorCode ierr,(*f)(PC,PCCompositeType);
341: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);
342: if (f) {
343: (*f)(pc,type);
344: }
345: return(0);
346: }
350: /*@
351: PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
352: for alphaI + R + S
353:
354: Collective on PC
356: Input Parameter:
357: + pc - the preconditioner context
358: - alpha - scale on identity
360: Level: Developer
362: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
363: @*/
364: PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
365: {
366: PetscErrorCode ierr,(*f)(PC,PetscScalar);
370: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);
371: if (f) {
372: (*f)(pc,alpha);
373: }
374: return(0);
375: }
379: /*@C
380: PCCompositeAddPC - Adds another PC to the composite PC.
381:
382: Collective on PC
384: Input Parameters:
385: + pc - the preconditioner context
386: - type - the type of the new preconditioner
388: Level: Developer
390: .keywords: PC, composite preconditioner, add
391: @*/
392: PetscErrorCode PCCompositeAddPC(PC pc,PCType type)
393: {
394: PetscErrorCode ierr,(*f)(PC,PCType);
398: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);
399: if (f) {
400: (*f)(pc,type);
401: }
402: return(0);
403: }
407: /*@
408: PCCompositeGetPC - Gets one of the PC objects in the composite PC.
409:
410: Not Collective
412: Input Parameter:
413: + pc - the preconditioner context
414: - n - the number of the pc requested
416: Output Parameters:
417: . subpc - the PC requested
419: Level: Developer
421: .keywords: PC, get, composite preconditioner, sub preconditioner
423: .seealso: PCCompositeAddPC()
424: @*/
425: PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
426: {
427: PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
432: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);
433: if (f) {
434: (*f)(pc,n,subpc);
435: } else {
436: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
437: }
438: return(0);
439: }
443: /*@
444: PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
445: the matrix used to define the preconditioner) is used to compute
446: the residual when the multiplicative scheme is used.
448: Collective on PC
450: Input Parameters:
451: . pc - the preconditioner context
453: Options Database Key:
454: . -pc_composite_true - Activates PCCompositeSetUseTrue()
456: Note:
457: For the common case in which the preconditioning and linear
458: system matrices are identical, this routine is unnecessary.
460: Level: Developer
462: .keywords: PC, composite preconditioner, set, true, flag
464: .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
465: @*/
466: PetscErrorCode PCCompositeSetUseTrue(PC pc)
467: {
468: PetscErrorCode ierr,(*f)(PC);
472: PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);
473: if (f) {
474: (*f)(pc);
475: }
476: return(0);
477: }
479: /* -------------------------------------------------------------------------------------------*/
481: /*MC
482: PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
484: Options Database Keys:
485: . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
486: . -pc_composite_true - Activates PCCompositeSetUseTrue()
488: Level: intermediate
490: Concepts: composing solvers
492: Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
493: inner PCs to be PCKSP.
494: Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
495: the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
498: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
499: PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
500: PCCompositeGetPC(), PCCompositeSetUseTrue()
502: M*/
507: PetscErrorCode PCCreate_Composite(PC pc)
508: {
510: PC_Composite *jac;
513: PetscNew(PC_Composite,&jac);
514: PetscLogObjectMemory(pc,sizeof(PC_Composite));
515: pc->ops->apply = PCApply_Composite_Additive;
516: pc->ops->setup = PCSetUp_Composite;
517: pc->ops->destroy = PCDestroy_Composite;
518: pc->ops->setfromoptions = PCSetFromOptions_Composite;
519: pc->ops->view = PCView_Composite;
520: pc->ops->applyrichardson = 0;
522: pc->data = (void*)jac;
523: jac->type = PC_COMPOSITE_ADDITIVE;
524: jac->work1 = 0;
525: jac->work2 = 0;
526: jac->head = 0;
528: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
529: PCCompositeSetType_Composite);
530: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
531: PCCompositeAddPC_Composite);
532: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
533: PCCompositeGetPC_Composite);
534: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
535: PCCompositeSetUseTrue_Composite);
536: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
537: PCCompositeSpecialSetAlpha_Composite);
539: return(0);
540: }