Actual source code: rgpolygon.c
slepc-3.6.1 2015-09-03
1: /*
2: Region defined by a set of vertices.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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: #define VERTMAX 30
28: typedef struct {
29: PetscInt n; /* number of vertices */
30: PetscScalar *vr,*vi; /* array of vertices (vi not used in complex scalars) */
31: } RG_POLYGON;
35: static PetscErrorCode RGPolygonSetVertices_Polygon(RG rg,PetscInt n,PetscScalar *vr,PetscScalar *vi)
36: {
38: PetscInt i;
39: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
42: if (n<3) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"At least 3 vertices required, you provided %s",n);
43: if (n>VERTMAX) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Too many points, maximum allowed is %d",VERTMAX);
44: if (ctx->n) {
45: PetscFree(ctx->vr);
46: #if !defined(PETSC_USE_COMPLEX)
47: PetscFree(ctx->vi);
48: #endif
49: }
50: ctx->n = n;
51: PetscMalloc1(n,&ctx->vr);
52: #if !defined(PETSC_USE_COMPLEX)
53: PetscMalloc1(n,&ctx->vi);
54: #endif
55: for (i=0;i<n;i++) {
56: ctx->vr[i] = vr[i];
57: #if !defined(PETSC_USE_COMPLEX)
58: ctx->vi[i] = vi[i];
59: #endif
60: }
61: return(0);
62: }
66: /*@
67: RGPolygonSetVertices - Sets the vertices that define the polygon region.
69: Logically Collective on RG
71: Input Parameters:
72: + rg - the region context
73: . n - number of vertices
74: . vr - array of vertices
75: - vi - array of vertices (imaginary part)
77: Options Database Keys:
78: + -rg_polygon_vertices - Sets the vertices
79: - -rg_polygon_verticesi - Sets the vertices (imaginary part)
81: Notes:
82: In the case of complex scalars, only argument vr is used, containing
83: the complex vertices; the list of vertices can be provided in the
84: command line with a comma-separated list of complex values
85: [+/-][realnumber][+/-]realnumberi with no spaces.
87: When PETSc is built with real scalars, the real and imaginary parts of
88: the vertices must be provided in two separate arrays (or two lists in
89: the command line).
91: Level: advanced
93: .seealso: RGPolygonGetVertices()
94: @*/
95: PetscErrorCode RGPolygonSetVertices(RG rg,PetscInt n,PetscScalar *vr,PetscScalar *vi)
96: {
103: #if !defined(PETSC_USE_COMPLEX)
105: #endif
106: PetscTryMethod(rg,"RGPolygonSetVertices_C",(RG,PetscInt,PetscScalar*,PetscScalar*),(rg,n,vr,vi));
107: return(0);
108: }
112: static PetscErrorCode RGPolygonGetVertices_Polygon(RG rg,PetscInt *n,PetscScalar **vr,PetscScalar **vi)
113: {
114: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
117: if (n) *n = ctx->n;
118: if (vr) *vr = ctx->vr;
119: if (vi) *vi = ctx->vi;
120: return(0);
121: }
125: /*@
126: RGPolygonGetVertices - Gets the vertices that define the polygon region.
128: Not Collective
130: Input Parameter:
131: . rg - the region context
133: Output Parameters:
134: + n - number of vertices
135: . vr - array of vertices
136: - vi - array of vertices (imaginary part)
138: Notes:
139: The returned arrays must NOT be freed by the calling application.
141: Level: advanced
143: .seealso: RGPolygonSetVertices()
144: @*/
145: PetscErrorCode RGPolygonGetVertices(RG rg,PetscInt *n,PetscScalar **vr,PetscScalar **vi)
146: {
151: PetscTryMethod(rg,"RGPolygonGetVertices_C",(RG,PetscInt*,PetscScalar**,PetscScalar**),(rg,n,vr,vi));
152: return(0);
153: }
157: PetscErrorCode RGView_Polygon(RG rg,PetscViewer viewer)
158: {
160: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
161: PetscBool isascii;
162: PetscInt i;
163: char str[50];
166: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
167: if (isascii) {
168: PetscViewerASCIIPrintf(viewer,"vertices: ");
169: for (i=0;i<ctx->n;i++) {
170: #if defined(PETSC_USE_COMPLEX)
171: SlepcSNPrintfScalar(str,50,ctx->vr[i],PETSC_FALSE);
172: #else
173: if (ctx->vi[i]!=0.0) {
174: PetscSNPrintf(str,50,"%g%+gi",(double)ctx->vr[i],(double)ctx->vi[i]);
175: } else {
176: PetscSNPrintf(str,50,"%g",(double)ctx->vr[i]);
177: }
178: #endif
179: PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->n-1)?",":"");
180: }
181: PetscViewerASCIIPrintf(viewer,"\n");
182: }
183: return(0);
184: }
188: PetscErrorCode RGIsTrivial_Polygon(RG rg,PetscBool *trivial)
189: {
190: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
193: *trivial = (ctx->n)? PETSC_FALSE: PETSC_TRUE;
194: return(0);
195: }
199: PetscErrorCode RGComputeContour_Polygon(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
200: {
201: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
204: if (!ctx->n) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONGSTATE,"No vertices have been set yet");
205: SETERRQ(PetscObjectComm((PetscObject)rg),1,"Not implemented yet");
206: return(0);
207: }
211: PetscErrorCode RGCheckInside_Polygon(RG rg,PetscReal px,PetscReal py,PetscInt *inout)
212: {
213: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
214: PetscReal val,x[VERTMAX],y[VERTMAX];
215: PetscBool mx,my,nx,ny;
216: PetscInt i,j;
219: for (i=0;i<ctx->n;i++) {
220: #if defined(PETSC_USE_COMPLEX)
221: x[i] = PetscRealPart(ctx->vr[i])-px;
222: y[i] = PetscImaginaryPart(ctx->vr[i])-py;
223: #else
224: x[i] = ctx->vr[i]-px;
225: y[i] = ctx->vi[i]-py;
226: #endif
227: }
228: *inout = -1;
229: for (i=0;i<ctx->n;i++) {
230: j = (i+1)%ctx->n;
231: mx = (x[i]>=0.0)? PETSC_TRUE: PETSC_FALSE;
232: nx = (x[j]>=0.0)? PETSC_TRUE: PETSC_FALSE;
233: my = (y[i]>=0.0)? PETSC_TRUE: PETSC_FALSE;
234: ny = (y[j]>=0.0)? PETSC_TRUE: PETSC_FALSE;
235: if (!((my||ny) && (mx||nx)) || (mx&&nx)) continue;
236: if (((my && ny && (mx||nx)) && (!(mx&&nx)))) {
237: *inout = -*inout;
238: continue;
239: }
240: val = (y[i]*x[j]-x[i]*y[j])/(x[j]-x[i]);
241: if (PetscAbs(val)<10*PETSC_MACHINE_EPSILON) {
242: *inout = 0;
243: return(0);
244: } else if (val>0.0) *inout = -*inout;
245: }
246: return(0);
247: }
251: PetscErrorCode RGSetFromOptions_Polygon(PetscOptions *PetscOptionsObject,RG rg)
252: {
254: PetscScalar array[VERTMAX];
255: PetscInt i,k;
256: PetscBool flg,flgi=PETSC_FALSE;
257: #if !defined(PETSC_USE_COMPLEX)
258: PetscScalar arrayi[VERTMAX];
259: PetscInt ki;
260: #else
261: PetscScalar *arrayi=NULL;
262: #endif
265: PetscOptionsHead(PetscOptionsObject,"RG Polygon Options");
267: k = VERTMAX;
268: for (i=0;i<k;i++) array[i] = 0;
269: PetscOptionsScalarArray("-rg_polygon_vertices","Vertices of polygon","RGPolygonSetVertices",array,&k,&flg);
270: #if !defined(PETSC_USE_COMPLEX)
271: ki = VERTMAX;
272: for (i=0;i<ki;i++) arrayi[i] = 0;
273: PetscOptionsScalarArray("-rg_polygon_verticesi","Vertices of polygon (imaginary part)","RGPolygonSetVertices",arrayi,&ki,&flgi);
274: 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);
275: #endif
276: if (flg || flgi) {
277: RGPolygonSetVertices(rg,k,array,arrayi);
278: }
280: PetscOptionsTail();
281: return(0);
282: }
286: PetscErrorCode RGDestroy_Polygon(RG rg)
287: {
289: RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
292: if (ctx->n) {
293: PetscFree(ctx->vr);
294: #if !defined(PETSC_USE_COMPLEX)
295: PetscFree(ctx->vi);
296: #endif
297: }
298: PetscFree(rg->data);
299: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",NULL);
300: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",NULL);
301: return(0);
302: }
306: PETSC_EXTERN PetscErrorCode RGCreate_Polygon(RG rg)
307: {
308: RG_POLYGON *polygon;
312: PetscNewLog(rg,&polygon);
313: rg->data = (void*)polygon;
315: rg->ops->istrivial = RGIsTrivial_Polygon;
316: rg->ops->computecontour = RGComputeContour_Polygon;
317: rg->ops->checkinside = RGCheckInside_Polygon;
318: rg->ops->setfromoptions = RGSetFromOptions_Polygon;
319: rg->ops->view = RGView_Polygon;
320: rg->ops->destroy = RGDestroy_Polygon;
321: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",RGPolygonSetVertices_Polygon);
322: PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",RGPolygonGetVertices_Polygon);
323: return(0);
324: }