Actual source code: galerkin.c
2: /*
3: Defines a preconditioner defined by R^T S R
4: */
5: #include <private/pcimpl.h> /*I "petscpc.h" I*/
6: #include <petscksp.h> /*I "petscksp.h" I*/
8: typedef struct {
9: KSP ksp;
10: Mat R,P;
11: Vec b,x;
12: } PC_Galerkin;
16: static PetscErrorCode PCApply_Galerkin(PC pc,Vec x,Vec y)
17: {
18: PetscErrorCode ierr;
19: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
22: if (jac->R) {MatRestrict(jac->R,x,jac->b);}
23: else {MatRestrict(jac->P,x,jac->b);}
24: KSPSolve(jac->ksp,jac->b,jac->x);
25: if (jac->P) {MatInterpolate(jac->P,jac->x,y);}
26: else {MatInterpolate(jac->R,jac->x,y);}
27: return(0);
28: }
32: static PetscErrorCode PCSetUp_Galerkin(PC pc)
33: {
34: PetscErrorCode ierr;
35: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
36: PetscBool a;
37: Vec *xx,*yy;
40: if (!jac->x) {
41: KSPGetOperatorsSet(jac->ksp,&a,PETSC_NULL);
42: if (!a) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set operator of PCGALERKIN KSP with PCGalerkinGetKSP()/KSPSetOperators()");
43: KSPGetVecs(jac->ksp,1,&xx,1,&yy);
44: jac->x = *xx;
45: jac->b = *yy;
46: PetscFree(xx);
47: PetscFree(yy);
48: }
49: if (!jac->R && !jac->P) {
50: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or interpolation of PCGALERKIN with PCGalerkinSetRestriction/Interpolation()");
51: }
52: /* should check here that sizes of R/P match size of a */
53: return(0);
54: }
58: static PetscErrorCode PCReset_Galerkin(PC pc)
59: {
60: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
61: PetscErrorCode ierr;
64: MatDestroy(&jac->R);
65: MatDestroy(&jac->P);
66: VecDestroy(&jac->x);
67: VecDestroy(&jac->b);
68: KSPReset(jac->ksp);
69: return(0);
70: }
74: static PetscErrorCode PCDestroy_Galerkin(PC pc)
75: {
76: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
77: PetscErrorCode ierr;
80: PCReset_Galerkin(pc);
81: KSPDestroy(&jac->ksp);
82: PetscFree(pc->data);
83: return(0);
84: }
88: static PetscErrorCode PCView_Galerkin(PC pc,PetscViewer viewer)
89: {
90: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
91: PetscErrorCode ierr;
92: PetscBool iascii;
95: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
96: if (iascii) {
97: PetscViewerASCIIPrintf(viewer,"Galerkin PC\n");
98: PetscViewerASCIIPrintf(viewer,"KSP on Galerkin follow\n");
99: PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
100: } else {
101: SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGalerkin",((PetscObject)viewer)->type_name);
102: }
103: KSPView(jac->ksp,viewer);
104: return(0);
105: }
110: PetscErrorCode PCGalerkinGetKSP_Galerkin(PC pc,KSP *ksp)
111: {
112: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
115: *ksp = jac->ksp;
116: return(0);
117: }
123: PetscErrorCode PCGalerkinSetRestriction_Galerkin(PC pc,Mat R)
124: {
125: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
126: PetscErrorCode ierr;
129: PetscObjectReference((PetscObject)R);
130: MatDestroy(&jac->R);
131: jac->R = R;
132: return(0);
133: }
139: PetscErrorCode PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P)
140: {
141: PC_Galerkin *jac = (PC_Galerkin*)pc->data;
142: PetscErrorCode ierr;
145: PetscObjectReference((PetscObject)P);
146: MatDestroy(&jac->P);
147: jac->P = P;
148: return(0);
149: }
152: /* -------------------------------------------------------------------------------- */
155: /*@
156: PCGalerkinSetRestriction - Sets the restriction operator for the "Galerkin-type" preconditioner
157:
158: Logically Collective on PC
160: Input Parameter:
161: + pc - the preconditioner context
162: - R - the restriction operator
164: Notes: Either this or PCGalerkinSetInterpolation() or both must be called
166: Level: Intermediate
168: .keywords: PC, set, Galerkin preconditioner
170: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
171: PCGalerkinSetInterpolation(), PCGalerkinGetKSP()
173: @*/
174: PetscErrorCode PCGalerkinSetRestriction(PC pc,Mat R)
175: {
180: PetscTryMethod(pc,"PCGalerkinSetRestriction_C",(PC,Mat),(pc,R));
181: return(0);
182: }
186: /*@
187: PCGalerkinSetInterpolation - Sets the interpolation operator for the "Galerkin-type" preconditioner
188:
189: Logically Collective on PC
191: Input Parameter:
192: + pc - the preconditioner context
193: - R - the interpolation operator
195: Notes: Either this or PCGalerkinSetRestriction() or both must be called
197: Level: Intermediate
199: .keywords: PC, set, Galerkin preconditioner
201: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
202: PCGalerkinSetRestriction(), PCGalerkinGetKSP()
204: @*/
205: PetscErrorCode PCGalerkinSetInterpolation(PC pc,Mat P)
206: {
211: PetscTryMethod(pc,"PCGalerkinSetInterpolation_C",(PC,Mat),(pc,P));
212: return(0);
213: }
217: /*@
218: PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC.
219:
220: Not Collective
222: Input Parameter:
223: . pc - the preconditioner context
225: Output Parameters:
226: . ksp - the KSP object
228: Level: Intermediate
230: .keywords: PC, get, Galerkin preconditioner, sub preconditioner
232: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN,
233: PCGalerkinSetRestriction(), PCGalerkinSetInterpolation()
235: @*/
236: PetscErrorCode PCGalerkinGetKSP(PC pc,KSP *ksp)
237: {
243: PetscUseMethod(pc,"PCGalerkinGetKSP_C",(PC,KSP *),(pc,ksp));
244: return(0);
245: }
248: /* -------------------------------------------------------------------------------------------*/
250: /*MC
251: PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T)
253: $ Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by
254: $ PCGalerkinGetKSP(pc,&ksp); KSPSetOperations(ksp,A,....)
256: Level: intermediate
258: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
259: PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP()
261: M*/
266: PetscErrorCode PCCreate_Galerkin(PC pc)
267: {
269: PC_Galerkin *jac;
272: PetscNewLog(pc,PC_Galerkin,&jac);
273: pc->ops->apply = PCApply_Galerkin;
274: pc->ops->setup = PCSetUp_Galerkin;
275: pc->ops->reset = PCReset_Galerkin;
276: pc->ops->destroy = PCDestroy_Galerkin;
277: pc->ops->view = PCView_Galerkin;
278: pc->ops->applyrichardson = 0;
280: KSPCreate(((PetscObject)pc)->comm,&jac->ksp);
281: PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1);
283: pc->data = (void*)jac;
285: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinSetRestriction_C","PCGalerkinSetRestriction_Galerkin",
286: PCGalerkinSetRestriction_Galerkin);
287: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinSetInterpolation_C","PCGalerkinSetInterpolation_Galerkin",
288: PCGalerkinSetInterpolation_Galerkin);
289: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinGetKSP_C","PCGalerkinGetKSP_Galerkin",
290: PCGalerkinGetKSP_Galerkin);
291: return(0);
292: }