Actual source code: rgring.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    Ring region, similar to the ellipse but with a start and end angle,
  3:    together with the width.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <slepc/private/rgimpl.h>      /*I "slepcrg.h" I*/

 27: typedef struct {
 28:   PetscScalar center;     /* center of the ellipse */
 29:   PetscReal   radius;     /* radius of the ellipse */
 30:   PetscReal   vscale;     /* vertical scale of the ellipse */
 31:   PetscReal   start_ang;  /* start angle */
 32:   PetscReal   end_ang;    /* end angle */
 33:   PetscReal   width;      /* ring width */
 34: } RG_RING;

 38: static PetscErrorCode RGRingSetParameters_Ring(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale,PetscReal start_ang,PetscReal end_ang,PetscReal width)
 39: {
 40:   RG_RING *ctx = (RG_RING*)rg->data;

 43:   ctx->center = center;
 44:   if (radius == PETSC_DEFAULT) {
 45:     ctx->radius = 1.0;
 46:   } else {
 47:     if (radius<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The radius argument must be > 0.0");
 48:     ctx->radius = radius;
 49:   }
 50:   if (vscale<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The vscale argument must be > 0.0");
 51:   ctx->vscale = vscale;
 52:   if (start_ang == PETSC_DEFAULT) {
 53:     ctx->start_ang = 0.0;
 54:   } else {
 55:     if (start_ang<0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The right-hand side angle argument must be >= 0.0");
 56:     if (start_ang>1.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The right-hand side angle argument must be <= 1.0");
 57:     ctx->start_ang = start_ang;
 58:   }
 59:   if (end_ang == PETSC_DEFAULT) {
 60:     ctx->end_ang = 1.0;
 61:   } else {
 62:     if (end_ang<0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The left-hand side angle argument must be >= 0.0");
 63:     if (end_ang>1.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The left-hand side angle argument must be <= 1.0");
 64:     ctx->end_ang = end_ang;
 65:   }
 66:   if (ctx->start_ang>ctx->end_ang) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The right-hand side angle argument must be smaller than left one");
 67:   if (width == PETSC_DEFAULT) {
 68:     ctx->width = 0.1;
 69:   } else {
 70:     if (width<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The width argument must be > 0.0");
 71:     ctx->width = width;
 72:   }
 73:   return(0);
 74: }

 78: /*@
 79:    RGRingSetParameters - Sets the parameters defining the ring region.

 81:    Logically Collective on RG

 83:    Input Parameters:
 84: +  rg        - the region context
 85: .  center    - center of the ellipse
 86: .  radius    - radius of the ellipse
 87: .  vscale    - vertical scale of the ellipse
 88: .  start_ang - the right-hand side angle
 89: .  end_ang   - the left-hand side angle
 90: -  width     - width of the ring

 92:    Options Database Keys:
 93: +  -rg_ring_center     - Sets the center
 94: .  -rg_ring_radius     - Sets the radius
 95: .  -rg_ring_vscale     - Sets the vertical scale
 96: .  -rg_ring_startangle - Sets the right-hand side angle
 97: .  -rg_ring_endangle   - Sets the left-hand side angle
 98: -  -rg_ring_width      - Sets the width of the ring

100:    Notes:
101:    The values of center, radius and vscale have the same meaning as in the
102:    ellipse region. The startangle and endangle define the span of the ring
103:    (by default it is the whole ring), while the width is the separation
104:    between the two concentric ellipses (above and below the radius by
105:    width/2). The start and end angles are expressed as a fraction of the
106:    circumference: the allowed range is [0..1], with 0 corresponding to 0
107:    radians, 0.25 to pi/2 radians, and so on.

109:    In the case of complex scalars, a complex center can be provided in the
110:    command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
111:    -rg_ring_center 1.0+2.0i

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

115:    Level: advanced

117: .seealso: RGRingGetParameters()
118: @*/
119: PetscErrorCode RGRingSetParameters(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale,PetscReal start_ang,PetscReal end_ang,PetscReal width)
120: {

131:   PetscTryMethod(rg,"RGRingSetParameters_C",(RG,PetscScalar,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal),(rg,center,radius,vscale,start_ang,end_ang,width));
132:   return(0);
133: }

137: static PetscErrorCode RGRingGetParameters_Ring(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale,PetscReal *start_ang,PetscReal *end_ang,PetscReal *width)
138: {
139:   RG_RING *ctx = (RG_RING*)rg->data;

142:   if (center)    *center    = ctx->center;
143:   if (radius)    *radius    = ctx->radius;
144:   if (vscale)    *vscale    = ctx->vscale;
145:   if (start_ang) *start_ang = ctx->start_ang;
146:   if (end_ang)   *end_ang   = ctx->end_ang;
147:   if (width)     *width     = ctx->width;
148:   return(0);
149: }

153: /*@
154:    RGRingGetParameters - Gets the parameters that define the ring region.

156:    Not Collective

158:    Input Parameter:
159: .  rg     - the region context

161:    Output Parameters:
162: +  center    - center of the region
163: .  radius    - radius of the region
164: .  vscale    - vertical scale of the region
165: .  start_ang - the right-hand side angle
166: .  end_ang   - the left-hand side angle
167: -  width     - width of the ring

169:    Level: advanced

171: .seealso: RGRingSetParameters()
172: @*/
173: PetscErrorCode RGRingGetParameters(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale,PetscReal *start_ang,PetscReal *end_ang,PetscReal *width)
174: {

179:   PetscTryMethod(rg,"RGRingGetParameters_C",(RG,PetscScalar*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(rg,center,radius,vscale,start_ang,end_ang,width));
180:   return(0);
181: }

185: PetscErrorCode RGView_Ring(RG rg,PetscViewer viewer)
186: {
188:   RG_RING        *ctx = (RG_RING*)rg->data;
189:   PetscBool      isascii;
190:   char           str[50];

193:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
194:   if (isascii) {
195:     SlepcSNPrintfScalar(str,50,ctx->center,PETSC_FALSE);
196:     PetscViewerASCIIPrintf(viewer,"center: %s, radius: %g, vscale: %g, start angle: %g, end angle: %g, ring width: %g\n",str,RGShowReal(ctx->radius),RGShowReal(ctx->vscale),ctx->start_ang,ctx->end_ang,ctx->width);
197:   }
198:   return(0);
199: }

203: PetscErrorCode RGIsTrivial_Ring(RG rg,PetscBool *trivial)
204: {
205:   RG_RING *ctx = (RG_RING*)rg->data;

208:   if (rg->complement) *trivial = (ctx->radius==0.0)? PETSC_TRUE: PETSC_FALSE;
209:   else *trivial = (ctx->radius>=PETSC_MAX_REAL)? PETSC_TRUE: PETSC_FALSE;
210:   return(0);
211: }

215: PetscErrorCode RGComputeContour_Ring(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
216: {
218:   SETERRQ(PetscObjectComm((PetscObject)rg),1,"Not implemented yet");
219:   return(0);
220: }

224: PetscErrorCode RGCheckInside_Ring(RG rg,PetscReal px,PetscReal py,PetscInt *inside)
225: {
226:   RG_RING   *ctx = (RG_RING*)rg->data;
227:   PetscReal dx,dy,r;

230:   /* outer ellipse */
231: #if defined(PETSC_USE_COMPLEX)
232:   dx = (px-PetscRealPart(ctx->center))/(ctx->radius+ctx->width/2.0);
233:   dy = (py-PetscImaginaryPart(ctx->center))/(ctx->radius+ctx->width/2.0);
234: #else
235:   dx = (px-ctx->center)/(ctx->radius+ctx->width/2.0);
236:   dy = py/(ctx->radius+ctx->width/2.0);
237: #endif
238:   r = 1.0-dx*dx-(dy*dy)/(ctx->vscale*ctx->vscale);
239:   *inside = PetscSign(r);
240:   /* inner ellipse */
241: #if defined(PETSC_USE_COMPLEX)
242:   dx = (px-PetscRealPart(ctx->center))/(ctx->radius-ctx->width/2.0);
243:   dy = (py-PetscImaginaryPart(ctx->center))/(ctx->radius-ctx->width/2.0);
244: #else
245:   dx = (px-ctx->center)/(ctx->radius-ctx->width/2.0);
246:   dy = py/(ctx->radius-ctx->width/2.0);
247: #endif
248:   r = -1.0+dx*dx+(dy*dy)/(ctx->vscale*ctx->vscale);
249:   *inside *= PetscSign(r);
250:   /* check angles */
251: #if defined(PETSC_USE_COMPLEX)
252:   dx = (px-PetscRealPart(ctx->center));
253:   dy = (py-PetscImaginaryPart(ctx->center));
254: #else
255:   dx = px-ctx->center;
256:   dy = py;
257: #endif
258:   if (dx == 0) {
259:     if (dy == 0) r = -1;
260:     else if (dy > 0) r = 0.25;
261:     else r = 0.75;
262:   } else if (dx > 0) {
263:     r = PetscAtanReal((dy/ctx->vscale)/dx);
264:     if (dy >= 0) r /= 2*PETSC_PI;
265:     else r = r/(2*PETSC_PI)+1;
266:   } else r = PetscAtanReal((dy/ctx->vscale)/dx)/(2*PETSC_PI)+0.5;
267:   if (r>=ctx->start_ang && r<=ctx->end_ang && *inside == 1) *inside = 1;
268:   else *inside = 0;
269:   return(0);
270: }

274: PetscErrorCode RGSetFromOptions_Ring(PetscOptions *PetscOptionsObject,RG rg)
275: {
277:   PetscScalar    s;
278:   PetscReal      r1,r2,r3,r4,r5;
279:   PetscBool      flg1,flg2,flg3,flg4,flg5,flg6;

282:   PetscOptionsHead(PetscOptionsObject,"RG Ring Options");

284:   RGRingGetParameters(rg,&s,&r1,&r2,&r3,&r4,&r5);
285:   PetscOptionsScalar("-rg_ring_center","Center of ellipse","RGRingSetParameters",s,&s,&flg1);
286:   PetscOptionsReal("-rg_ring_radius","Radius of ellipse","RGRingSetParameters",r1,&r1,&flg2);
287:   PetscOptionsReal("-rg_ring_vscale","Vertical scale of ellipse","RGRingSetParameters",r2,&r2,&flg3);
288:   PetscOptionsReal("-rg_ring_startangle","Right-hand side angle","RGRingSetParameters",r3,&r3,&flg4);
289:   PetscOptionsReal("-rg_ring_endangle","Left-hand side angle","RGRingSetParameters",r4,&r4,&flg5);
290:   PetscOptionsReal("-rg_ring_width","Width of ring","RGRingSetParameters",r5,&r5,&flg6);
291:   if (flg1 || flg2 || flg3 || flg4 || flg5 || flg6) {
292:     RGRingSetParameters(rg,s,r1,r2,r3,r4,r5);
293:   }

295:   PetscOptionsTail();
296:   return(0);
297: }

301: PetscErrorCode RGDestroy_Ring(RG rg)
302: {

306:   PetscFree(rg->data);
307:   PetscObjectComposeFunction((PetscObject)rg,"RGRingSetParameters_C",NULL);
308:   PetscObjectComposeFunction((PetscObject)rg,"RGRingGetParameters_C",NULL);
309:   return(0);
310: }

314: PETSC_EXTERN PetscErrorCode RGCreate_Ring(RG rg)
315: {
316:   RG_RING        *ring;

320:   PetscNewLog(rg,&ring);
321:   ring->center    = 0.0;
322:   ring->radius    = 1.0;
323:   ring->vscale    = 1.0;
324:   ring->start_ang = 0.0;
325:   ring->end_ang   = 1.0;
326:   ring->width     = 0.1;
327:   rg->data = (void*)ring;

329:   rg->ops->istrivial      = RGIsTrivial_Ring;
330:   rg->ops->computecontour = RGComputeContour_Ring;
331:   rg->ops->checkinside    = RGCheckInside_Ring;
332:   rg->ops->setfromoptions = RGSetFromOptions_Ring;
333:   rg->ops->view           = RGView_Ring;
334:   rg->ops->destroy        = RGDestroy_Ring;
335:   PetscObjectComposeFunction((PetscObject)rg,"RGRingSetParameters_C",RGRingSetParameters_Ring);
336:   PetscObjectComposeFunction((PetscObject)rg,"RGRingGetParameters_C",RGRingGetParameters_Ring);
337:   return(0);
338: }