Actual source code: rgellipse.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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>
 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;

 27:   ctx->center = center;
 28:   if (radius == PETSC_DEFAULT) {
 29:     ctx->radius = 1.0;
 30:   } else {
 32:     ctx->radius = radius;
 33:   }
 35:   ctx->vscale = vscale;
 36:   PetscFunctionReturn(0);
 37: }

 39: /*@
 40:    RGEllipseSetParameters - Sets the parameters defining the ellipse region.

 42:    Logically Collective on rg

 44:    Input Parameters:
 45: +  rg     - the region context
 46: .  center - center of the ellipse
 47: .  radius - radius of the ellipse
 48: -  vscale - vertical scale of the ellipse

 50:    Options Database Keys:
 51: +  -rg_ellipse_center - Sets the center
 52: .  -rg_ellipse_radius - Sets the radius
 53: -  -rg_ellipse_vscale - Sets the vertical scale

 55:    Notes:
 56:    In the case of complex scalars, a complex center can be provided in the
 57:    command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
 58:    -rg_ellipse_center 1.0+2.0i

 60:    When PETSc is built with real scalars, the center is restricted to a real value.

 62:    Level: advanced

 64: .seealso: RGEllipseGetParameters()
 65: @*/
 66: PetscErrorCode RGEllipseSetParameters(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale)
 67: {
 72:   PetscTryMethod(rg,"RGEllipseSetParameters_C",(RG,PetscScalar,PetscReal,PetscReal),(rg,center,radius,vscale));
 73:   PetscFunctionReturn(0);
 74: }

 76: static PetscErrorCode RGEllipseGetParameters_Ellipse(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
 77: {
 78:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;

 80:   if (center) *center = ctx->center;
 81:   if (radius) *radius = ctx->radius;
 82:   if (vscale) *vscale = ctx->vscale;
 83:   PetscFunctionReturn(0);
 84: }

 86: /*@C
 87:    RGEllipseGetParameters - Gets the parameters that define the ellipse region.

 89:    Not Collective

 91:    Input Parameter:
 92: .  rg     - the region context

 94:    Output Parameters:
 95: +  center - center of the region
 96: .  radius - radius of the region
 97: -  vscale - vertical scale of the region

 99:    Level: advanced

101: .seealso: RGEllipseSetParameters()
102: @*/
103: PetscErrorCode RGEllipseGetParameters(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale)
104: {
106:   PetscUseMethod(rg,"RGEllipseGetParameters_C",(RG,PetscScalar*,PetscReal*,PetscReal*),(rg,center,radius,vscale));
107:   PetscFunctionReturn(0);
108: }

110: PetscErrorCode RGView_Ellipse(RG rg,PetscViewer viewer)
111: {
112:   RG_ELLIPSE     *ctx = (RG_ELLIPSE*)rg->data;
113:   PetscBool      isdraw,isascii;
114:   int            winw,winh;
115:   PetscDraw      draw;
116:   PetscDrawAxis  axis;
117:   PetscReal      cx,cy,r,ab,cd,lx,ly,w,scale=1.2;
118:   char           str[50];

120:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
121:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
122:   if (isascii) {
123:     SlepcSNPrintfScalar(str,sizeof(str),ctx->center,PETSC_FALSE);
124:     PetscViewerASCIIPrintf(viewer,"  center: %s, radius: %g, vscale: %g\n",str,RGShowReal(ctx->radius),RGShowReal(ctx->vscale));
125:   } else if (isdraw) {
126:     PetscViewerDrawGetDraw(viewer,0,&draw);
127:     PetscDrawCheckResizedWindow(draw);
128:     PetscDrawGetWindowSize(draw,&winw,&winh);
129:     winw = PetscMax(winw,1); winh = PetscMax(winh,1);
130:     PetscDrawClear(draw);
131:     PetscDrawSetTitle(draw,"Ellipse region");
132:     PetscDrawAxisCreate(draw,&axis);
133:     cx = PetscRealPart(ctx->center)*rg->sfactor;
134:     cy = PetscImaginaryPart(ctx->center)*rg->sfactor;
135:     r  = ctx->radius*rg->sfactor;
136:     lx = 2*r;
137:     ly = 2*r*ctx->vscale;
138:     ab = cx;
139:     cd = cy;
140:     w  = scale*PetscMax(lx/winw,ly/winh)/2;
141:     PetscDrawAxisSetLimits(axis,ab-w*winw,ab+w*winw,cd-w*winh,cd+w*winh);
142:     PetscDrawAxisDraw(axis);
143:     PetscDrawAxisDestroy(&axis);
144:     PetscDrawEllipse(draw,cx,cy,2*r,2*r*ctx->vscale,PETSC_DRAW_RED);
145:     PetscDrawFlush(draw);
146:     PetscDrawSave(draw);
147:     PetscDrawPause(draw);
148:   }
149:   PetscFunctionReturn(0);
150: }

152: PetscErrorCode RGIsTrivial_Ellipse(RG rg,PetscBool *trivial)
153: {
154:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;

156:   if (rg->complement) *trivial = PetscNot(ctx->radius);
157:   else *trivial = PetscNot(ctx->radius<PETSC_MAX_REAL);
158:   PetscFunctionReturn(0);
159: }

161: PetscErrorCode RGComputeContour_Ellipse(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
162: {
163:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
164:   PetscReal  theta;
165:   PetscInt   i;

167:   for (i=0;i<n;i++) {
168:     theta = 2.0*PETSC_PI*(i+0.5)/n;
169: #if defined(PETSC_USE_COMPLEX)
170:     cr[i] = ctx->center + ctx->radius*PetscCMPLX(PetscCosReal(theta),ctx->vscale*PetscSinReal(theta));
171: #else
172:     if (cr) cr[i] = ctx->center + ctx->radius*PetscCosReal(theta);
173:     if (ci) ci[i] = ctx->radius*ctx->vscale*PetscSinReal(theta);
174: #endif
175:   }
176:   PetscFunctionReturn(0);
177: }

179: PetscErrorCode RGComputeBoundingBox_Ellipse(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
180: {
181:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;

183:   if (a) *a = PetscRealPart(ctx->center) - ctx->radius;
184:   if (b) *b = PetscRealPart(ctx->center) + ctx->radius;
185:   if (c) *c = PetscImaginaryPart(ctx->center) - ctx->radius*ctx->vscale;
186:   if (d) *d = PetscImaginaryPart(ctx->center) + ctx->radius*ctx->vscale;
187:   PetscFunctionReturn(0);
188: }

190: PetscErrorCode RGComputeQuadrature_Ellipse(RG rg,RGQuadRule quad,PetscInt n,PetscScalar *z,PetscScalar *zn,PetscScalar *w)
191: {
192:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
193:   PetscReal  theta;
194:   PetscInt   i;

196:   for (i=0;i<n;i++) {
197: #if defined(PETSC_USE_COMPLEX)
198:     theta = 2.0*PETSC_PI*(i+0.5)/n;
199:     zn[i] = PetscCMPLX(PetscCosReal(theta),ctx->vscale*PetscSinReal(theta));
200:     w[i]  = rg->sfactor*ctx->radius*(PetscCMPLX(ctx->vscale*PetscCosReal(theta),PetscSinReal(theta)))/n;
201: #else
202:     theta = PETSC_PI*(i+0.5)/n;
203:     zn[i] = PetscCosReal(theta);
204:     w[i]  = PetscCosReal((n-1)*theta)/n;
205:     z[i]  = rg->sfactor*(ctx->center + ctx->radius*zn[i]);
206: #endif
207:   }
208:   PetscFunctionReturn(0);
209: }

211: PetscErrorCode RGCheckInside_Ellipse(RG rg,PetscReal px,PetscReal py,PetscInt *inside)
212: {
213:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;
214:   PetscReal  dx,dy,r;

216: #if defined(PETSC_USE_COMPLEX)
217:   dx = (px-PetscRealPart(ctx->center))/ctx->radius;
218:   dy = (py-PetscImaginaryPart(ctx->center))/ctx->radius;
219: #else
220:   dx = (px-ctx->center)/ctx->radius;
221:   dy = py/ctx->radius;
222: #endif
223:   r = 1.0-dx*dx-(dy*dy)/(ctx->vscale*ctx->vscale);
224:   *inside = PetscSign(r);
225:   PetscFunctionReturn(0);
226: }

228: PetscErrorCode RGIsAxisymmetric_Ellipse(RG rg,PetscBool vertical,PetscBool *symm)
229: {
230:   RG_ELLIPSE *ctx = (RG_ELLIPSE*)rg->data;

232:   if (vertical) *symm = (PetscRealPart(ctx->center) == 0.0)? PETSC_TRUE: PETSC_FALSE;
233:   else *symm = (PetscImaginaryPart(ctx->center) == 0.0)? PETSC_TRUE: PETSC_FALSE;
234:   PetscFunctionReturn(0);
235: }

237: PetscErrorCode RGSetFromOptions_Ellipse(PetscOptionItems *PetscOptionsObject,RG rg)
238: {
239:   PetscScalar    s;
240:   PetscReal      r1,r2;
241:   PetscBool      flg1,flg2,flg3;

243:   PetscOptionsHead(PetscOptionsObject,"RG Ellipse Options");

245:     RGEllipseGetParameters(rg,&s,&r1,&r2);
246:     PetscOptionsScalar("-rg_ellipse_center","Center of ellipse","RGEllipseSetParameters",s,&s,&flg1);
247:     PetscOptionsReal("-rg_ellipse_radius","Radius of ellipse","RGEllipseSetParameters",r1,&r1,&flg2);
248:     PetscOptionsReal("-rg_ellipse_vscale","Vertical scale of ellipse","RGEllipseSetParameters",r2,&r2,&flg3);
249:     if (flg1 || flg2 || flg3) RGEllipseSetParameters(rg,s,r1,r2);

251:   PetscOptionsTail();
252:   PetscFunctionReturn(0);
253: }

255: PetscErrorCode RGDestroy_Ellipse(RG rg)
256: {
257:   PetscFree(rg->data);
258:   PetscObjectComposeFunction((PetscObject)rg,"RGEllipseSetParameters_C",NULL);
259:   PetscObjectComposeFunction((PetscObject)rg,"RGEllipseGetParameters_C",NULL);
260:   PetscFunctionReturn(0);
261: }

263: SLEPC_EXTERN PetscErrorCode RGCreate_Ellipse(RG rg)
264: {
265:   RG_ELLIPSE     *ellipse;

267:   PetscNewLog(rg,&ellipse);
268:   ellipse->center = 0.0;
269:   ellipse->radius = PETSC_MAX_REAL;
270:   ellipse->vscale = 1.0;
271:   rg->data = (void*)ellipse;

273:   rg->ops->istrivial         = RGIsTrivial_Ellipse;
274:   rg->ops->computecontour    = RGComputeContour_Ellipse;
275:   rg->ops->computebbox       = RGComputeBoundingBox_Ellipse;
276:   rg->ops->computequadrature = RGComputeQuadrature_Ellipse;
277:   rg->ops->checkinside       = RGCheckInside_Ellipse;
278:   rg->ops->isaxisymmetric    = RGIsAxisymmetric_Ellipse;
279:   rg->ops->setfromoptions    = RGSetFromOptions_Ellipse;
280:   rg->ops->view              = RGView_Ellipse;
281:   rg->ops->destroy           = RGDestroy_Ellipse;
282:   PetscObjectComposeFunction((PetscObject)rg,"RGEllipseSetParameters_C",RGEllipseSetParameters_Ellipse);
283:   PetscObjectComposeFunction((PetscObject)rg,"RGEllipseGetParameters_C",RGEllipseGetParameters_Ellipse);
284:   PetscFunctionReturn(0);
285: }