Actual source code: rgpolygon.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: Polygonal region defined by a set of vertices
12: */
14: #include <slepc/private/rgimpl.h> /*I "slepcrg.h" I*/
15: #include <petscdraw.h>
17: #define VERTMAX 30
19: typedef struct {
20: PetscInt n; /* number of vertices */
21: PetscScalar *vr,*vi; /* array of vertices (vi not used in complex scalars) */
22: } RG_POLYGON;
24: PetscErrorCode RGComputeBoundingBox_Polygon(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
26: #if !defined(PETSC_USE_COMPLEX)
27: static PetscBool CheckSymmetry(PetscInt n,PetscScalar *vr,PetscScalar *vi)
28: {
29: PetscInt i,j,k;
30: /* find change of sign in imaginary part */
31: j = vi[0]!=0.0? 0: 1;
32: for (k=j+1;k<n;k++) {
33: if (vi[k]!=0.0) {
34: if (vi[k]*vi[j]<0.0) break;
35: j++;
36: }
37: }
38: if (k==n) return (j==1)? PETSC_TRUE: PETSC_FALSE;
39: /* check pairing vertices */
40: for (i=0;i<n/2;i++) {
41: if (vr[k]!=vr[j] || vi[k]!=-vi[j]) return PETSC_FALSE;
42: k = (k+1)%n;
43: j = (j-1+n)%n;
44: }
45: return PETSC_TRUE;
46: }
47: #endif
49: static PetscErrorCode RGPolygonSetVertices_Polygon(RG rg,PetscInt n,PetscScalar *vr,PetscScalar *vi)
50: {
52: PetscInt i;
53: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
56: if (n<3) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"At least 3 vertices required, you provided %s",n);
57: if (n>VERTMAX) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Too many points, maximum allowed is %d",VERTMAX);
58: #if !defined(PETSC_USE_COMPLEX)
59: if (!CheckSymmetry(n,vr,vi)) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"In real scalars the region must be symmetric wrt real axis");
60: #endif
61: if (ctx->n) {
62: PetscFree(ctx->vr);
63: #if !defined(PETSC_USE_COMPLEX)
64: PetscFree(ctx->vi);
65: #endif
66: }
67: ctx->n = n;
68: PetscMalloc1(n,&ctx->vr);
69: #if !defined(PETSC_USE_COMPLEX)
70: PetscMalloc1(n,&ctx->vi);
71: #endif
72: for (i=0;i<n;i++) {
73: ctx->vr[i] = vr[i];
74: #if !defined(PETSC_USE_COMPLEX)
75: ctx->vi[i] = vi[i];
76: #endif
77: }
78: return(0);
79: }
81: /*@
82: RGPolygonSetVertices - Sets the vertices that define the polygon region.
84: Logically Collective on RG
86: Input Parameters:
87: + rg - the region context
88: . n - number of vertices
89: . vr - array of vertices
90: - vi - array of vertices (imaginary part)
92: Options Database Keys:
93: + -rg_polygon_vertices - Sets the vertices
94: - -rg_polygon_verticesi - Sets the vertices (imaginary part)
96: Notes:
97: In the case of complex scalars, only argument vr is used, containing
98: the complex vertices; the list of vertices can be provided in the
99: command line with a comma-separated list of complex values
100: [+/-][realnumber][+/-]realnumberi with no spaces.
102: When PETSc is built with real scalars, the real and imaginary parts of
103: the vertices must be provided in two separate arrays (or two lists in
104: the command line). In this case, the region must be symmetric with
105: respect to the real axis.
107: Level: advanced
109: .seealso: RGPolygonGetVertices()
110: @*/
111: PetscErrorCode RGPolygonSetVertices(RG rg,PetscInt n,PetscScalar vr[],PetscScalar vi[])
112: {
119: #if !defined(PETSC_USE_COMPLEX)
121: #endif
122: PetscTryMethod(rg,"RGPolygonSetVertices_C",(RG,PetscInt,PetscScalar*,PetscScalar*),(rg,n,vr,vi));
123: return(0);
124: }
126: static PetscErrorCode RGPolygonGetVertices_Polygon(RG rg,PetscInt *n,PetscScalar **vr,PetscScalar **vi)
127: {
128: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
131: if (n) *n = ctx->n;
132: if (vr) *vr = ctx->vr;
133: if (vi) *vi = ctx->vi;
134: return(0);
135: }
137: /*@
138: RGPolygonGetVertices - Gets the vertices that define the polygon region.
140: Not Collective
142: Input Parameter:
143: . rg - the region context
145: Output Parameters:
146: + n - number of vertices
147: . vr - array of vertices
148: - vi - array of vertices (imaginary part)
150: Notes:
151: The returned arrays must NOT be freed by the calling application.
153: Level: advanced
155: .seealso: RGPolygonSetVertices()
156: @*/
157: PetscErrorCode RGPolygonGetVertices(RG rg,PetscInt *n,PetscScalar **vr,PetscScalar **vi)
158: {
163: PetscUseMethod(rg,"RGPolygonGetVertices_C",(RG,PetscInt*,PetscScalar**,PetscScalar**),(rg,n,vr,vi));
164: return(0);
165: }
167: PetscErrorCode RGView_Polygon(RG rg,PetscViewer viewer)
168: {
170: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
171: PetscBool isdraw,isascii;
172: int winw,winh;
173: PetscDraw draw;
174: PetscDrawAxis axis;
175: PetscReal a,b,c,d,ab,cd,lx,ly,w,x0,y0,x1,y1,scale=1.2;
176: PetscInt i;
177: char str[50];
180: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
181: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
182: if (isascii) {
183: PetscViewerASCIIPrintf(viewer," vertices: ");
184: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
185: for (i=0;i<ctx->n;i++) {
186: #if defined(PETSC_USE_COMPLEX)
187: SlepcSNPrintfScalar(str,50,ctx->vr[i],PETSC_FALSE);
188: #else
189: if (ctx->vi[i]!=0.0) {
190: PetscSNPrintf(str,50,"%g%+gi",(double)ctx->vr[i],(double)ctx->vi[i]);
191: } else {
192: PetscSNPrintf(str,50,"%g",(double)ctx->vr[i]);
193: }
194: #endif
195: PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->n-1)?", ":"");
196: }
197: PetscViewerASCIIPrintf(viewer,"\n");
198: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
199: } else if (isdraw) {
200: PetscViewerDrawGetDraw(viewer,0,&draw);
201: PetscDrawCheckResizedWindow(draw);
202: PetscDrawGetWindowSize(draw,&winw,&winh);
203: winw = PetscMax(winw,1); winh = PetscMax(winh,1);
204: PetscDrawClear(draw);
205: PetscDrawSetTitle(draw,"Polygonal region");
206: PetscDrawAxisCreate(draw,&axis);
207: RGComputeBoundingBox_Polygon(rg,&a,&b,&c,&d);
208: a *= rg->sfactor;
209: b *= rg->sfactor;
210: c *= rg->sfactor;
211: d *= rg->sfactor;
212: lx = b-a;
213: ly = d-c;
214: ab = (a+b)/2;
215: cd = (c+d)/2;
216: w = scale*PetscMax(lx/winw,ly/winh)/2;
217: PetscDrawAxisSetLimits(axis,ab-w*winw,ab+w*winw,cd-w*winh,cd+w*winh);
218: PetscDrawAxisDraw(axis);
219: PetscDrawAxisDestroy(&axis);
220: for (i=0;i<ctx->n;i++) {
221: #if defined(PETSC_USE_COMPLEX)
222: x0 = PetscRealPart(ctx->vr[i]); y0 = PetscImaginaryPart(ctx->vr[i]);
223: if (i<ctx->n-1) {
224: x1 = PetscRealPart(ctx->vr[i+1]); y1 = PetscImaginaryPart(ctx->vr[i+1]);
225: } else {
226: x1 = PetscRealPart(ctx->vr[0]); y1 = PetscImaginaryPart(ctx->vr[0]);
227: }
228: #else
229: x0 = ctx->vr[i]; y0 = ctx->vi[i];
230: if (i<ctx->n-1) {
231: x1 = ctx->vr[i+1]; y1 = ctx->vi[i+1];
232: } else {
233: x1 = ctx->vr[0]; y1 = ctx->vi[0];
234: }
235: #endif
236: PetscDrawLine(draw,x0*rg->sfactor,y0*rg->sfactor,x1*rg->sfactor,y1*rg->sfactor,PETSC_DRAW_MAGENTA);
237: }
238: PetscDrawFlush(draw);
239: PetscDrawSave(draw);
240: PetscDrawPause(draw);
241: }
242: return(0);
243: }
245: PetscErrorCode RGIsTrivial_Polygon(RG rg,PetscBool *trivial)
246: {
247: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
250: *trivial = PetscNot(ctx->n);
251: return(0);
252: }
254: PetscErrorCode RGComputeContour_Polygon(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
255: {
256: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
257: PetscReal length,h,d,rem=0.0;
258: PetscInt k=1,idx=ctx->n-1,i;
259: PetscBool ini=PETSC_FALSE;
260: PetscScalar incr;
261: #if !defined(PETSC_USE_COMPLEX)
262: PetscScalar inci;
263: #endif
266: if (!ctx->n) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONGSTATE,"No vertices have been set yet");
267: length = SlepcAbsEigenvalue(ctx->vr[0]-ctx->vr[ctx->n-1],ctx->vi[0]-ctx->vi[ctx->n-1]);
268: for (i=0;i<ctx->n-1;i++) length += SlepcAbsEigenvalue(ctx->vr[i]-ctx->vr[i+1],ctx->vi[i]-ctx->vi[i+1]);
269: h = length/n;
270: cr[0] = ctx->vr[0];
271: #if !defined(PETSC_USE_COMPLEX)
272: ci[0] = ctx->vi[0];
273: #endif
274: incr = ctx->vr[ctx->n-1]-ctx->vr[0];
275: #if !defined(PETSC_USE_COMPLEX)
276: inci = ctx->vi[ctx->n-1]-ctx->vi[0];
277: #endif
278: d = SlepcAbsEigenvalue(incr,inci);
279: incr /= d;
280: #if !defined(PETSC_USE_COMPLEX)
281: inci /= d;
282: #endif
283: while (k<n) {
284: if (ini) {
285: incr = ctx->vr[idx]-ctx->vr[idx+1];
286: #if !defined(PETSC_USE_COMPLEX)
287: inci = ctx->vi[idx]-ctx->vi[idx+1];
288: #endif
289: d = SlepcAbsEigenvalue(incr,inci);
290: incr /= d;
291: #if !defined(PETSC_USE_COMPLEX)
292: inci /= d;
293: #endif
294: if (rem+d>h) {
295: cr[k] = ctx->vr[idx+1]+incr*(h-rem);
296: #if !defined(PETSC_USE_COMPLEX)
297: ci[k] = ctx->vi[idx+1]+inci*(h-rem);
298: #endif
299: k++;
300: ini = PETSC_FALSE;
301: } else {rem += d; idx--;}
302: } else {
303: #if !defined(PETSC_USE_COMPLEX)
304: rem = SlepcAbsEigenvalue(ctx->vr[idx]-cr[k-1],ctx->vi[idx]-ci[k-1]);
305: #else
306: rem = PetscAbsScalar(ctx->vr[idx]-cr[k-1]);
307: #endif
308: if (rem>h) {
309: cr[k] = cr[k-1]+incr*h;
310: #if !defined(PETSC_USE_COMPLEX)
311: ci[k] = ci[k-1]+inci*h;
312: #endif
313: k++;
314: } else {ini = PETSC_TRUE; idx--;}
315: }
316: }
317: return(0);
318: }
320: PetscErrorCode RGComputeBoundingBox_Polygon(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
321: {
322: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
323: PetscInt i;
326: *a = PETSC_MAX_REAL;
327: *b = -PETSC_MAX_REAL;
328: *c = PETSC_MAX_REAL;
329: *d = -PETSC_MAX_REAL;
330: for (i=0;i<ctx->n;i++) {
331: #if defined(PETSC_USE_COMPLEX)
332: if (a) *a = PetscMin(*a,PetscRealPart(ctx->vr[i]));
333: if (b) *b = PetscMax(*b,PetscRealPart(ctx->vr[i]));
334: if (c) *c = PetscMin(*c,PetscImaginaryPart(ctx->vr[i]));
335: if (d) *d = PetscMax(*d,PetscImaginaryPart(ctx->vr[i]));
336: #else
337: if (a) *a = PetscMin(*a,ctx->vr[i]);
338: if (b) *b = PetscMax(*b,ctx->vr[i]);
339: if (c) *c = PetscMin(*c,ctx->vi[i]);
340: if (d) *d = PetscMax(*d,ctx->vi[i]);
341: #endif
342: }
343: return(0);
344: }
346: PetscErrorCode RGCheckInside_Polygon(RG rg,PetscReal px,PetscReal py,PetscInt *inout)
347: {
348: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
349: PetscReal val,x[VERTMAX],y[VERTMAX];
350: PetscBool mx,my,nx,ny;
351: PetscInt i,j;
354: for (i=0;i<ctx->n;i++) {
355: #if defined(PETSC_USE_COMPLEX)
356: x[i] = PetscRealPart(ctx->vr[i])-px;
357: y[i] = PetscImaginaryPart(ctx->vr[i])-py;
358: #else
359: x[i] = ctx->vr[i]-px;
360: y[i] = ctx->vi[i]-py;
361: #endif
362: }
363: *inout = -1;
364: for (i=0;i<ctx->n;i++) {
365: j = (i+1)%ctx->n;
366: mx = PetscNot(x[i]<0.0);
367: nx = PetscNot(x[j]<0.0);
368: my = PetscNot(y[i]<0.0);
369: ny = PetscNot(y[j]<0.0);
370: if (!((my||ny) && (mx||nx)) || (mx&&nx)) continue;
371: if (((my && ny && (mx||nx)) && (!(mx&&nx)))) {
372: *inout = -*inout;
373: continue;
374: }
375: val = (y[i]*x[j]-x[i]*y[j])/(x[j]-x[i]);
376: if (PetscAbs(val)<10*PETSC_MACHINE_EPSILON) {
377: *inout = 0;
378: return(0);
379: } else if (val>0.0) *inout = -*inout;
380: }
381: return(0);
382: }
384: PetscErrorCode RGSetFromOptions_Polygon(PetscOptionItems *PetscOptionsObject,RG rg)
385: {
387: PetscScalar array[VERTMAX];
388: PetscInt i,k;
389: PetscBool flg,flgi=PETSC_FALSE;
390: #if !defined(PETSC_USE_COMPLEX)
391: PetscScalar arrayi[VERTMAX];
392: PetscInt ki;
393: #else
394: PetscScalar *arrayi=NULL;
395: #endif
398: PetscOptionsHead(PetscOptionsObject,"RG Polygon Options");
400: k = VERTMAX;
401: for (i=0;i<k;i++) array[i] = 0;
402: PetscOptionsScalarArray("-rg_polygon_vertices","Vertices of polygon","RGPolygonSetVertices",array,&k,&flg);
403: #if !defined(PETSC_USE_COMPLEX)
404: ki = VERTMAX;
405: for (i=0;i<ki;i++) arrayi[i] = 0;
406: PetscOptionsScalarArray("-rg_polygon_verticesi","Vertices of polygon (imaginary part)","RGPolygonSetVertices",arrayi,&ki,&flgi);
407: 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);
408: #endif
409: if (flg || flgi) { RGPolygonSetVertices(rg,k,array,arrayi); }
411: PetscOptionsTail();
412: return(0);
413: }
415: PetscErrorCode RGDestroy_Polygon(RG rg)
416: {
418: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
421: if (ctx->n) {
422: PetscFree(ctx->vr);
423: #if !defined(PETSC_USE_COMPLEX)
424: PetscFree(ctx->vi);
425: #endif
426: }
427: PetscFree(rg->data);
428: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",NULL);
429: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",NULL);
430: return(0);
431: }
433: SLEPC_EXTERN PetscErrorCode RGCreate_Polygon(RG rg)
434: {
435: RG_POLYGON *polygon;
439: PetscNewLog(rg,&polygon);
440: rg->data = (void*)polygon;
442: rg->ops->istrivial = RGIsTrivial_Polygon;
443: rg->ops->computecontour = RGComputeContour_Polygon;
444: rg->ops->computebbox = RGComputeBoundingBox_Polygon;
445: rg->ops->checkinside = RGCheckInside_Polygon;
446: rg->ops->setfromoptions = RGSetFromOptions_Polygon;
447: rg->ops->view = RGView_Polygon;
448: rg->ops->destroy = RGDestroy_Polygon;
449: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",RGPolygonSetVertices_Polygon);
450: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",RGPolygonGetVertices_Polygon);
451: return(0);
452: }