Actual source code: rgpolygon.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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:    Polygonal region defined by a set of vertices
 12: */

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

 16: #define VERTMAX 30

 18: typedef struct {
 19:   PetscInt    n;         /* number of vertices */
 20:   PetscScalar *vr,*vi;   /* array of vertices (vi not used in complex scalars) */
 21: } RG_POLYGON;

 23: #if !defined(PETSC_USE_COMPLEX)
 24: static PetscBool CheckSymmetry(PetscInt n,PetscScalar *vr,PetscScalar *vi)
 25: {
 26:   PetscInt i,j,k;
 27:   /* find change of sign in imaginary part */
 28:   j = vi[0]!=0.0? 0: 1;
 29:   for (k=j+1;k<n;k++) {
 30:     if (vi[k]!=0.0) {
 31:       if (vi[k]*vi[j]<0.0) break;
 32:       j++;
 33:     }
 34:   }
 35:   if (k==n) return (j==1)? PETSC_TRUE: PETSC_FALSE;
 36:   /* check pairing vertices */
 37:   for (i=0;i<n/2;i++) {
 38:     if (vr[k]!=vr[j] || vi[k]!=-vi[j]) return PETSC_FALSE;
 39:     k = (k+1)%n;
 40:     j = (j-1+n)%n;
 41:   }
 42:   return PETSC_TRUE;
 43: }
 44: #endif

 46: static PetscErrorCode RGPolygonSetVertices_Polygon(RG rg,PetscInt n,PetscScalar *vr,PetscScalar *vi)
 47: {
 49:   PetscInt       i;
 50:   RG_POLYGON     *ctx = (RG_POLYGON*)rg->data;

 53:   if (n<3) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"At least 3 vertices required, you provided %s",n);
 54:   if (n>VERTMAX) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Too many points, maximum allowed is %d",VERTMAX);
 55: #if !defined(PETSC_USE_COMPLEX)
 56:   if (!CheckSymmetry(n,vr,vi)) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"In real scalars the region must be symmetric wrt real axis");
 57: #endif
 58:   if (ctx->n) {
 59:     PetscFree(ctx->vr);
 60: #if !defined(PETSC_USE_COMPLEX)
 61:     PetscFree(ctx->vi);
 62: #endif
 63:   }
 64:   ctx->n = n;
 65:   PetscMalloc1(n,&ctx->vr);
 66: #if !defined(PETSC_USE_COMPLEX)
 67:   PetscMalloc1(n,&ctx->vi);
 68: #endif
 69:   for (i=0;i<n;i++) {
 70:     ctx->vr[i] = vr[i];
 71: #if !defined(PETSC_USE_COMPLEX)
 72:     ctx->vi[i] = vi[i];
 73: #endif
 74:   }
 75:   return(0);
 76: }

 78: /*@
 79:    RGPolygonSetVertices - Sets the vertices that define the polygon region.

 81:    Logically Collective on RG

 83:    Input Parameters:
 84: +  rg - the region context
 85: .  n  - number of vertices
 86: .  vr - array of vertices
 87: -  vi - array of vertices (imaginary part)

 89:    Options Database Keys:
 90: +  -rg_polygon_vertices - Sets the vertices
 91: -  -rg_polygon_verticesi - Sets the vertices (imaginary part)

 93:    Notes:
 94:    In the case of complex scalars, only argument vr is used, containing
 95:    the complex vertices; the list of vertices can be provided in the
 96:    command line with a comma-separated list of complex values
 97:    [+/-][realnumber][+/-]realnumberi with no spaces.

 99:    When PETSc is built with real scalars, the real and imaginary parts of
100:    the vertices must be provided in two separate arrays (or two lists in
101:    the command line). In this case, the region must be symmetric with
102:    respect to the real axis.

104:    Level: advanced

106: .seealso: RGPolygonGetVertices()
107: @*/
108: PetscErrorCode RGPolygonSetVertices(RG rg,PetscInt n,PetscScalar vr[],PetscScalar vi[])
109: {

116: #if !defined(PETSC_USE_COMPLEX)
118: #endif
119:   PetscTryMethod(rg,"RGPolygonSetVertices_C",(RG,PetscInt,PetscScalar*,PetscScalar*),(rg,n,vr,vi));
120:   return(0);
121: }

123: static PetscErrorCode RGPolygonGetVertices_Polygon(RG rg,PetscInt *n,PetscScalar **vr,PetscScalar **vi)
124: {
125:   RG_POLYGON *ctx = (RG_POLYGON*)rg->data;

128:   if (n)  *n  = ctx->n;
129:   if (vr) *vr = ctx->vr;
130:   if (vi) *vi = ctx->vi;
131:   return(0);
132: }

134: /*@
135:    RGPolygonGetVertices - Gets the vertices that define the polygon region.

137:    Not Collective

139:    Input Parameter:
140: .  rg     - the region context

142:    Output Parameters:
143: +  n  - number of vertices
144: .  vr - array of vertices
145: -  vi - array of vertices (imaginary part)

147:    Notes:
148:    The returned arrays must NOT be freed by the calling application.

150:    Level: advanced

152: .seealso: RGPolygonSetVertices()
153: @*/
154: PetscErrorCode RGPolygonGetVertices(RG rg,PetscInt *n,PetscScalar **vr,PetscScalar **vi)
155: {

160:   PetscUseMethod(rg,"RGPolygonGetVertices_C",(RG,PetscInt*,PetscScalar**,PetscScalar**),(rg,n,vr,vi));
161:   return(0);
162: }

164: PetscErrorCode RGView_Polygon(RG rg,PetscViewer viewer)
165: {
167:   RG_POLYGON     *ctx = (RG_POLYGON*)rg->data;
168:   PetscBool      isascii;
169:   PetscInt       i;
170:   char           str[50];

173:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
174:   if (isascii) {
175:     PetscViewerASCIIPrintf(viewer,"  vertices: ");
176:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
177:     for (i=0;i<ctx->n;i++) {
178: #if defined(PETSC_USE_COMPLEX)
179:       SlepcSNPrintfScalar(str,50,ctx->vr[i],PETSC_FALSE);
180: #else
181:       if (ctx->vi[i]!=0.0) {
182:         PetscSNPrintf(str,50,"%g%+gi",(double)ctx->vr[i],(double)ctx->vi[i]);
183:       } else {
184:         PetscSNPrintf(str,50,"%g",(double)ctx->vr[i]);
185:       }
186: #endif
187:       PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->n-1)?", ":"");
188:     }
189:     PetscViewerASCIIPrintf(viewer,"\n");
190:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
191:   }
192:   return(0);
193: }

195: PetscErrorCode RGIsTrivial_Polygon(RG rg,PetscBool *trivial)
196: {
197:   RG_POLYGON *ctx = (RG_POLYGON*)rg->data;

200:   *trivial = PetscNot(ctx->n);
201:   return(0);
202: }

204: PetscErrorCode RGComputeContour_Polygon(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
205: {
206:   RG_POLYGON  *ctx = (RG_POLYGON*)rg->data;
207:   PetscReal   length,h,d,rem=0.0;
208:   PetscInt    k=1,idx=ctx->n-1,i;
209:   PetscBool   ini=PETSC_FALSE;
210:   PetscScalar incr;
211: #if !defined(PETSC_USE_COMPLEX)
212:   PetscScalar inci;
213: #endif

216:   if (!ctx->n) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONGSTATE,"No vertices have been set yet");
217:   length = SlepcAbsEigenvalue(ctx->vr[0]-ctx->vr[ctx->n-1],ctx->vi[0]-ctx->vi[ctx->n-1]);
218:   for (i=0;i<ctx->n-1;i++) length += SlepcAbsEigenvalue(ctx->vr[i]-ctx->vr[i+1],ctx->vi[i]-ctx->vi[i+1]);
219:   h = length/n;
220:   cr[0] = ctx->vr[0];
221: #if !defined(PETSC_USE_COMPLEX)
222:   ci[0] = ctx->vi[0];
223: #endif
224:   incr = ctx->vr[ctx->n-1]-ctx->vr[0];
225: #if !defined(PETSC_USE_COMPLEX)
226:   inci = ctx->vi[ctx->n-1]-ctx->vi[0];
227: #endif
228:   d = SlepcAbsEigenvalue(incr,inci);
229:   incr /= d;
230: #if !defined(PETSC_USE_COMPLEX)
231:   inci /= d;
232: #endif
233:   while (k<n) {
234:     if (ini) {
235:       incr = ctx->vr[idx]-ctx->vr[idx+1];
236: #if !defined(PETSC_USE_COMPLEX)
237:       inci = ctx->vi[idx]-ctx->vi[idx+1];
238: #endif
239:       d = SlepcAbsEigenvalue(incr,inci);
240:       incr /= d;
241: #if !defined(PETSC_USE_COMPLEX)
242:       inci /= d;
243: #endif
244:       if (rem+d>h) {
245:         cr[k] = ctx->vr[idx+1]+incr*(h-rem);
246: #if !defined(PETSC_USE_COMPLEX)
247:         ci[k] = ctx->vi[idx+1]+inci*(h-rem);
248: #endif
249:         k++;
250:         ini = PETSC_FALSE;
251:       } else {rem += d; idx--;}
252:     } else {
253: #if !defined(PETSC_USE_COMPLEX)
254:       rem = SlepcAbsEigenvalue(ctx->vr[idx]-cr[k-1],ctx->vi[idx]-ci[k-1]);
255: #else
256:       rem = PetscAbsScalar(ctx->vr[idx]-cr[k-1]);
257: #endif
258:       if (rem>h) {
259:         cr[k] = cr[k-1]+incr*h;
260: #if !defined(PETSC_USE_COMPLEX)
261:         ci[k] = ci[k-1]+inci*h;
262: #endif
263:         k++;
264:       } else {ini = PETSC_TRUE; idx--;}
265:     }
266:   }
267:   return(0);
268: }

270: PetscErrorCode RGComputeBoundingBox_Polygon(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
271: {
272:   RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
273:   PetscInt   i;

276:   *a =  PETSC_MAX_REAL;
277:   *b = -PETSC_MAX_REAL;
278:   *c =  PETSC_MAX_REAL;
279:   *d = -PETSC_MAX_REAL;
280:   for (i=0;i<ctx->n;i++) {
281: #if defined(PETSC_USE_COMPLEX)
282:     if (a) *a = PetscMin(*a,PetscRealPart(ctx->vr[i]));
283:     if (b) *b = PetscMax(*b,PetscRealPart(ctx->vr[i]));
284:     if (c) *c = PetscMin(*c,PetscImaginaryPart(ctx->vr[i]));
285:     if (d) *d = PetscMax(*d,PetscImaginaryPart(ctx->vr[i]));
286: #else
287:     if (a) *a = PetscMin(*a,ctx->vr[i]);
288:     if (b) *b = PetscMax(*b,ctx->vr[i]);
289:     if (c) *c = PetscMin(*c,ctx->vi[i]);
290:     if (d) *d = PetscMax(*d,ctx->vi[i]);
291: #endif
292:   }
293:   return(0);
294: }

296: PetscErrorCode RGCheckInside_Polygon(RG rg,PetscReal px,PetscReal py,PetscInt *inout)
297: {
298:   RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
299:   PetscReal  val,x[VERTMAX],y[VERTMAX];
300:   PetscBool  mx,my,nx,ny;
301:   PetscInt   i,j;

304:   for (i=0;i<ctx->n;i++) {
305: #if defined(PETSC_USE_COMPLEX)
306:     x[i] = PetscRealPart(ctx->vr[i])-px;
307:     y[i] = PetscImaginaryPart(ctx->vr[i])-py;
308: #else
309:     x[i] = ctx->vr[i]-px;
310:     y[i] = ctx->vi[i]-py;
311: #endif
312:   }
313:   *inout = -1;
314:   for (i=0;i<ctx->n;i++) {
315:     j = (i+1)%ctx->n;
316:     mx = PetscNot(x[i]<0.0);
317:     nx = PetscNot(x[j]<0.0);
318:     my = PetscNot(y[i]<0.0);
319:     ny = PetscNot(y[j]<0.0);
320:     if (!((my||ny) && (mx||nx)) || (mx&&nx)) continue;
321:     if (((my && ny && (mx||nx)) && (!(mx&&nx)))) {
322:       *inout = -*inout;
323:       continue;
324:     }
325:     val = (y[i]*x[j]-x[i]*y[j])/(x[j]-x[i]);
326:     if (PetscAbs(val)<10*PETSC_MACHINE_EPSILON) {
327:       *inout = 0;
328:       return(0);
329:     } else if (val>0.0) *inout = -*inout;
330:   }
331:   return(0);
332: }

334: PetscErrorCode RGSetFromOptions_Polygon(PetscOptionItems *PetscOptionsObject,RG rg)
335: {
337:   PetscScalar    array[VERTMAX];
338:   PetscInt       i,k;
339:   PetscBool      flg,flgi=PETSC_FALSE;
340: #if !defined(PETSC_USE_COMPLEX)
341:   PetscScalar    arrayi[VERTMAX];
342:   PetscInt       ki;
343: #else
344:   PetscScalar    *arrayi=NULL;
345: #endif

348:   PetscOptionsHead(PetscOptionsObject,"RG Polygon Options");

350:     k = VERTMAX;
351:     for (i=0;i<k;i++) array[i] = 0;
352:     PetscOptionsScalarArray("-rg_polygon_vertices","Vertices of polygon","RGPolygonSetVertices",array,&k,&flg);
353: #if !defined(PETSC_USE_COMPLEX)
354:     ki = VERTMAX;
355:     for (i=0;i<ki;i++) arrayi[i] = 0;
356:     PetscOptionsScalarArray("-rg_polygon_verticesi","Vertices of polygon (imaginary part)","RGPolygonSetVertices",arrayi,&ki,&flgi);
357:     if (ki!=k) SETERRQ2(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_SIZ,"The number of real %D and imaginary %D parts do not match",k,ki);
358: #endif
359:     if (flg || flgi) { RGPolygonSetVertices(rg,k,array,arrayi); }

361:   PetscOptionsTail();
362:   return(0);
363: }

365: PetscErrorCode RGDestroy_Polygon(RG rg)
366: {
368:   RG_POLYGON     *ctx = (RG_POLYGON*)rg->data;

371:   if (ctx->n) {
372:     PetscFree(ctx->vr);
373: #if !defined(PETSC_USE_COMPLEX)
374:     PetscFree(ctx->vi);
375: #endif
376:   }
377:   PetscFree(rg->data);
378:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",NULL);
379:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",NULL);
380:   return(0);
381: }

383: PETSC_EXTERN PetscErrorCode RGCreate_Polygon(RG rg)
384: {
385:   RG_POLYGON     *polygon;

389:   PetscNewLog(rg,&polygon);
390:   rg->data = (void*)polygon;

392:   rg->ops->istrivial      = RGIsTrivial_Polygon;
393:   rg->ops->computecontour = RGComputeContour_Polygon;
394:   rg->ops->computebbox    = RGComputeBoundingBox_Polygon;
395:   rg->ops->checkinside    = RGCheckInside_Polygon;
396:   rg->ops->setfromoptions = RGSetFromOptions_Polygon;
397:   rg->ops->view           = RGView_Polygon;
398:   rg->ops->destroy        = RGDestroy_Polygon;
399:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",RGPolygonSetVertices_Polygon);
400:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",RGPolygonGetVertices_Polygon);
401:   return(0);
402: }