Actual source code: rgellipse.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Region enclosed in an ellipse (aligned with the coordinate axes)
12: */
14: #include <slepc/private/rgimpl.h> /*I "slepcrg.h" I*/
15: #include <petscdraw.h>
17: typedef struct {
18: PetscScalar center; /* center of the ellipse */
19: PetscReal radius; /* radius of the ellipse */
20: PetscReal vscale; /* vertical scale of the ellipse */
21: } RG_ELLIPSE;
23: static PetscErrorCode RGEllipseSetParameters_Ellipse(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale)
24: {
25: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
28: ctx->center = center;
29: if (radius == PETSC_DEFAULT) {
30: ctx->radius = 1.0;
31: } else {
32: if (radius<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The radius argument must be > 0.0");
33: ctx->radius = radius;
34: }
35: if (vscale<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The vscale argument must be > 0.0");
36: ctx->vscale = vscale;
37: return(0);
38: }
40: /*@
41: RGEllipseSetParameters - Sets the parameters defining the ellipse region.
43: Logically Collective on RG
45: Input Parameters:
46: + rg - the region context
47: . center - center of the ellipse
48: . radius - radius of the ellipse
49: - vscale - vertical scale of the ellipse
51: Options Database Keys:
52: + -rg_ellipse_center - Sets the center
53: . -rg_ellipse_radius - Sets the radius
54: - -rg_ellipse_vscale - Sets the vertical scale
56: Notes:
57: In the case of complex scalars, a complex center can be provided in the
58: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
59: -rg_ellipse_center 1.0+2.0i
61: When PETSc is built with real scalars, the center is restricted to a real value.
63: Level: advanced
65: .seealso: RGEllipseGetParameters()
66: @*/
67: PetscErrorCode RGEllipseSetParameters(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale)
68: {
76: PetscTryMethod(rg,"RGEllipseSetParameters_C",(RG,PetscScalar,PetscReal,PetscReal),(rg,center,radius,vscale));
77: return(0);
78: }
80: static PetscErrorCode RGEllipseGetParameters_Ellipse(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
81: {
82: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
85: if (center) *center = ctx->center;
86: if (radius) *radius = ctx->radius;
87: if (vscale) *vscale = ctx->vscale;
88: return(0);
89: }
91: /*@
92: RGEllipseGetParameters - Gets the parameters that define the ellipse region.
94: Not Collective
96: Input Parameter:
97: . rg - the region context
99: Output Parameters:
100: + center - center of the region
101: . radius - radius of the region
102: - vscale - vertical scale of the region
104: Level: advanced
106: .seealso: RGEllipseSetParameters()
107: @*/
108: PetscErrorCode RGEllipseGetParameters(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
109: {
114: PetscUseMethod(rg,"RGEllipseGetParameters_C",(RG,PetscScalar*,PetscReal*,PetscReal*),(rg,center,radius,vscale));
115: return(0);
116: }
118: PetscErrorCode RGView_Ellipse(RG rg,PetscViewer viewer)
119: {
121: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
122: PetscBool isdraw,isascii;
123: int winw,winh;
124: PetscDraw draw;
125: PetscDrawAxis axis;
126: PetscReal cx,cy,r,ab,cd,lx,ly,w,scale=1.2;
127: char str[50];
130: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
131: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
132: if (isascii) {
133: SlepcSNPrintfScalar(str,50,ctx->center,PETSC_FALSE);
134: PetscViewerASCIIPrintf(viewer," center: %s, radius: %g, vscale: %g\n",str,RGShowReal(ctx->radius),RGShowReal(ctx->vscale));
135: } else if (isdraw) {
136: PetscViewerDrawGetDraw(viewer,0,&draw);
137: PetscDrawCheckResizedWindow(draw);
138: PetscDrawGetWindowSize(draw,&winw,&winh);
139: winw = PetscMax(winw,1); winh = PetscMax(winh,1);
140: PetscDrawClear(draw);
141: PetscDrawSetTitle(draw,"Ellipse region");
142: PetscDrawAxisCreate(draw,&axis);
143: cx = PetscRealPart(ctx->center)*rg->sfactor;
144: cy = PetscImaginaryPart(ctx->center)*rg->sfactor;
145: r = ctx->radius*rg->sfactor;
146: lx = 2*r;
147: ly = 2*r*ctx->vscale;
148: ab = cx;
149: cd = cy;
150: w = scale*PetscMax(lx/winw,ly/winh)/2;
151: PetscDrawAxisSetLimits(axis,ab-w*winw,ab+w*winw,cd-w*winh,cd+w*winh);
152: PetscDrawAxisDraw(axis);
153: PetscDrawAxisDestroy(&axis);
154: PetscDrawEllipse(draw,cx,cy,2*r,2*r*ctx->vscale,PETSC_DRAW_RED);
155: PetscDrawFlush(draw);
156: PetscDrawSave(draw);
157: PetscDrawPause(draw);
158: }
159: return(0);
160: }
162: PetscErrorCode RGIsTrivial_Ellipse(RG rg,PetscBool *trivial)
163: {
164: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
167: if (rg->complement) *trivial = PetscNot(ctx->radius);
168: else *trivial = PetscNot(ctx->radius<PETSC_MAX_REAL);
169: return(0);
170: }
172: PetscErrorCode RGComputeContour_Ellipse(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
173: {
174: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
175: PetscReal theta;
176: PetscInt i;
179: for (i=0;i<n;i++) {
180: theta = 2.0*PETSC_PI*(i+0.5)/n;
181: #if defined(PETSC_USE_COMPLEX)
182: cr[i] = ctx->center + ctx->radius*(PetscCosReal(theta)+ctx->vscale*PetscSinReal(theta)*PETSC_i);
183: #else
184: cr[i] = ctx->center + ctx->radius*PetscCosReal(theta);
185: ci[i] = ctx->radius*ctx->vscale*PetscSinReal(theta);
186: #endif
187: }
188: return(0);
189: }
191: PetscErrorCode RGComputeBoundingBox_Ellipse(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
192: {
193: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
196: if (a) *a = PetscRealPart(ctx->center) - ctx->radius;
197: if (b) *b = PetscRealPart(ctx->center) + ctx->radius;
198: if (c) *c = PetscImaginaryPart(ctx->center) - ctx->radius*ctx->vscale;
199: if (d) *d = PetscImaginaryPart(ctx->center) + ctx->radius*ctx->vscale;
200: return(0);
201: }
203: PetscErrorCode RGCheckInside_Ellipse(RG rg,PetscReal px,PetscReal py,PetscInt *inside)
204: {
205: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
206: PetscReal dx,dy,r;
209: #if defined(PETSC_USE_COMPLEX)
210: dx = (px-PetscRealPart(ctx->center))/ctx->radius;
211: dy = (py-PetscImaginaryPart(ctx->center))/ctx->radius;
212: #else
213: dx = (px-ctx->center)/ctx->radius;
214: dy = py/ctx->radius;
215: #endif
216: r = 1.0-dx*dx-(dy*dy)/(ctx->vscale*ctx->vscale);
217: *inside = PetscSign(r);
218: return(0);
219: }
221: PetscErrorCode RGSetFromOptions_Ellipse(PetscOptionItems *PetscOptionsObject,RG rg)
222: {
224: PetscScalar s;
225: PetscReal r1,r2;
226: PetscBool flg1,flg2,flg3;
229: PetscOptionsHead(PetscOptionsObject,"RG Ellipse Options");
231: RGEllipseGetParameters(rg,&s,&r1,&r2);
232: PetscOptionsScalar("-rg_ellipse_center","Center of ellipse","RGEllipseSetParameters",s,&s,&flg1);
233: PetscOptionsReal("-rg_ellipse_radius","Radius of ellipse","RGEllipseSetParameters",r1,&r1,&flg2);
234: PetscOptionsReal("-rg_ellipse_vscale","Vertical scale of ellipse","RGEllipseSetParameters",r2,&r2,&flg3);
235: if (flg1 || flg2 || flg3) { RGEllipseSetParameters(rg,s,r1,r2); }
237: PetscOptionsTail();
238: return(0);
239: }
241: PetscErrorCode RGDestroy_Ellipse(RG rg)
242: {
246: PetscFree(rg->data);
247: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseSetParameters_C",NULL);
248: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseGetParameters_C",NULL);
249: return(0);
250: }
252: SLEPC_EXTERN PetscErrorCode RGCreate_Ellipse(RG rg)
253: {
254: RG_ELLIPSE *ellipse;
258: PetscNewLog(rg,&ellipse);
259: ellipse->center = 0.0;
260: ellipse->radius = PETSC_MAX_REAL;
261: ellipse->vscale = 1.0;
262: rg->data = (void*)ellipse;
264: rg->ops->istrivial = RGIsTrivial_Ellipse;
265: rg->ops->computecontour = RGComputeContour_Ellipse;
266: rg->ops->computebbox = RGComputeBoundingBox_Ellipse;
267: rg->ops->checkinside = RGCheckInside_Ellipse;
268: rg->ops->setfromoptions = RGSetFromOptions_Ellipse;
269: rg->ops->view = RGView_Ellipse;
270: rg->ops->destroy = RGDestroy_Ellipse;
271: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseSetParameters_C",RGEllipseSetParameters_Ellipse);
272: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseGetParameters_C",RGEllipseGetParameters_Ellipse);
273: return(0);
274: }