Actual source code: rgellipse.c
slepc-3.7.0 2016-05-16
1: /*
2: Region enclosed in an ellipse (aligned with the coordinate axes).
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/rgimpl.h> /*I "slepcrg.h" I*/
26: typedef struct {
27: PetscScalar center; /* center of the ellipse */
28: PetscReal radius; /* radius of the ellipse */
29: PetscReal vscale; /* vertical scale of the ellipse */
30: } RG_ELLIPSE;
34: static PetscErrorCode RGEllipseSetParameters_Ellipse(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale)
35: {
36: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
39: ctx->center = center;
40: if (radius == PETSC_DEFAULT) {
41: ctx->radius = 1.0;
42: } else {
43: if (radius<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The radius argument must be > 0.0");
44: ctx->radius = radius;
45: }
46: if (vscale<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The vscale argument must be > 0.0");
47: ctx->vscale = vscale;
48: return(0);
49: }
53: /*@
54: RGEllipseSetParameters - Sets the parameters defining the ellipse region.
56: Logically Collective on RG
58: Input Parameters:
59: + rg - the region context
60: . center - center of the ellipse
61: . radius - radius of the ellipse
62: - vscale - vertical scale of the ellipse
64: Options Database Keys:
65: + -rg_ellipse_center - Sets the center
66: . -rg_ellipse_radius - Sets the radius
67: - -rg_ellipse_vscale - Sets the vertical scale
69: Notes:
70: In the case of complex scalars, a complex center can be provided in the
71: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
72: -rg_ellipse_center 1.0+2.0i
74: When PETSc is built with real scalars, the center is restricted to a real value.
76: Level: advanced
78: .seealso: RGEllipseGetParameters()
79: @*/
80: PetscErrorCode RGEllipseSetParameters(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale)
81: {
89: PetscTryMethod(rg,"RGEllipseSetParameters_C",(RG,PetscScalar,PetscReal,PetscReal),(rg,center,radius,vscale));
90: return(0);
91: }
95: static PetscErrorCode RGEllipseGetParameters_Ellipse(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
96: {
97: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
100: if (center) *center = ctx->center;
101: if (radius) *radius = ctx->radius;
102: if (vscale) *vscale = ctx->vscale;
103: return(0);
104: }
108: /*@
109: RGEllipseGetParameters - Gets the parameters that define the ellipse region.
111: Not Collective
113: Input Parameter:
114: . rg - the region context
116: Output Parameters:
117: + center - center of the region
118: . radius - radius of the region
119: - vscale - vertical scale of the region
121: Level: advanced
123: .seealso: RGEllipseSetParameters()
124: @*/
125: PetscErrorCode RGEllipseGetParameters(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
126: {
131: PetscUseMethod(rg,"RGEllipseGetParameters_C",(RG,PetscScalar*,PetscReal*,PetscReal*),(rg,center,radius,vscale));
132: return(0);
133: }
137: PetscErrorCode RGView_Ellipse(RG rg,PetscViewer viewer)
138: {
140: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
141: PetscBool isascii;
142: char str[50];
145: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
146: if (isascii) {
147: SlepcSNPrintfScalar(str,50,ctx->center,PETSC_FALSE);
148: PetscViewerASCIIPrintf(viewer,"center: %s, radius: %g, vscale: %g\n",str,RGShowReal(ctx->radius),RGShowReal(ctx->vscale));
149: }
150: return(0);
151: }
155: PetscErrorCode RGIsTrivial_Ellipse(RG rg,PetscBool *trivial)
156: {
157: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
160: if (rg->complement) *trivial = PetscNot(ctx->radius);
161: else *trivial = PetscNot(ctx->radius<PETSC_MAX_REAL);
162: return(0);
163: }
167: PetscErrorCode RGComputeContour_Ellipse(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
168: {
169: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
170: PetscReal theta;
171: PetscInt i;
174: for (i=0;i<n;i++) {
175: theta = 2.0*PETSC_PI*(i+0.5)/n;
176: #if defined(PETSC_USE_COMPLEX)
177: cr[i] = ctx->center + ctx->radius*(PetscCosReal(theta)+ctx->vscale*PetscSinReal(theta)*PETSC_i);
178: #else
179: cr[i] = ctx->center + ctx->radius*PetscCosReal(theta);
180: ci[i] = ctx->radius*ctx->vscale*PetscSinReal(theta);
181: #endif
182: }
183: return(0);
184: }
188: PetscErrorCode RGCheckInside_Ellipse(RG rg,PetscReal px,PetscReal py,PetscInt *inside)
189: {
190: RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
191: PetscReal dx,dy,r;
194: #if defined(PETSC_USE_COMPLEX)
195: dx = (px-PetscRealPart(ctx->center))/ctx->radius;
196: dy = (py-PetscImaginaryPart(ctx->center))/ctx->radius;
197: #else
198: dx = (px-ctx->center)/ctx->radius;
199: dy = py/ctx->radius;
200: #endif
201: r = 1.0-dx*dx-(dy*dy)/(ctx->vscale*ctx->vscale);
202: *inside = PetscSign(r);
203: return(0);
204: }
208: PetscErrorCode RGSetFromOptions_Ellipse(PetscOptionItems *PetscOptionsObject,RG rg)
209: {
211: PetscScalar s;
212: PetscReal r1,r2;
213: PetscBool flg1,flg2,flg3;
216: PetscOptionsHead(PetscOptionsObject,"RG Ellipse Options");
218: RGEllipseGetParameters(rg,&s,&r1,&r2);
219: PetscOptionsScalar("-rg_ellipse_center","Center of ellipse","RGEllipseSetParameters",s,&s,&flg1);
220: PetscOptionsReal("-rg_ellipse_radius","Radius of ellipse","RGEllipseSetParameters",r1,&r1,&flg2);
221: PetscOptionsReal("-rg_ellipse_vscale","Vertical scale of ellipse","RGEllipseSetParameters",r2,&r2,&flg3);
222: if (flg1 || flg2 || flg3) {
223: RGEllipseSetParameters(rg,s,r1,r2);
224: }
226: PetscOptionsTail();
227: return(0);
228: }
232: PetscErrorCode RGDestroy_Ellipse(RG rg)
233: {
237: PetscFree(rg->data);
238: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseSetParameters_C",NULL);
239: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseGetParameters_C",NULL);
240: return(0);
241: }
245: PETSC_EXTERN PetscErrorCode RGCreate_Ellipse(RG rg)
246: {
247: RG_ELLIPSE *ellipse;
251: PetscNewLog(rg,&ellipse);
252: ellipse->center = 0.0;
253: ellipse->radius = 1.0;
254: ellipse->vscale = 1.0;
255: rg->data = (void*)ellipse;
257: rg->ops->istrivial = RGIsTrivial_Ellipse;
258: rg->ops->computecontour = RGComputeContour_Ellipse;
259: rg->ops->checkinside = RGCheckInside_Ellipse;
260: rg->ops->setfromoptions = RGSetFromOptions_Ellipse;
261: rg->ops->view = RGView_Ellipse;
262: rg->ops->destroy = RGDestroy_Ellipse;
263: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseSetParameters_C",RGEllipseSetParameters_Ellipse);
264: PetscObjectComposeFunction((PetscObject)rg,"RGEllipseGetParameters_C",RGEllipseGetParameters_Ellipse);
265: return(0);
266: }