Actual source code: rgpolygon.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    Region defined by a set of vertices.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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:   PetscUseMethod(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:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
170:     for (i=0;i<ctx->n;i++) {
171: #if defined(PETSC_USE_COMPLEX)
172:       SlepcSNPrintfScalar(str,50,ctx->vr[i],PETSC_FALSE);
173: #else
174:       if (ctx->vi[i]!=0.0) {
175:         PetscSNPrintf(str,50,"%g%+gi",(double)ctx->vr[i],(double)ctx->vi[i]);
176:       } else {
177:         PetscSNPrintf(str,50,"%g",(double)ctx->vr[i]);
178:       }
179: #endif
180:       PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->n-1)?",":"");
181:     }
182:     PetscViewerASCIIPrintf(viewer,"\n");
183:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
184:   }
185:   return(0);
186: }

190: PetscErrorCode RGIsTrivial_Polygon(RG rg,PetscBool *trivial)
191: {
192:   RG_POLYGON *ctx = (RG_POLYGON*)rg->data;

195:   *trivial = PetscNot(ctx->n);
196:   return(0);
197: }

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

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

269: PetscErrorCode RGCheckInside_Polygon(RG rg,PetscReal px,PetscReal py,PetscInt *inout)
270: {
271:   RG_POLYGON *ctx = (RG_POLYGON*)rg->data;
272:   PetscReal  val,x[VERTMAX],y[VERTMAX];
273:   PetscBool  mx,my,nx,ny;
274:   PetscInt   i,j;

277:   for (i=0;i<ctx->n;i++) {
278: #if defined(PETSC_USE_COMPLEX)
279:     x[i] = PetscRealPart(ctx->vr[i])-px;
280:     y[i] = PetscImaginaryPart(ctx->vr[i])-py;
281: #else
282:     x[i] = ctx->vr[i]-px;
283:     y[i] = ctx->vi[i]-py;
284: #endif
285:   }
286:   *inout = -1;
287:   for (i=0;i<ctx->n;i++) {
288:     j = (i+1)%ctx->n;
289:     mx = PetscNot(x[i]<0.0);
290:     nx = PetscNot(x[j]<0.0);
291:     my = PetscNot(y[i]<0.0);
292:     ny = PetscNot(y[j]<0.0);
293:     if (!((my||ny) && (mx||nx)) || (mx&&nx)) continue;
294:     if (((my && ny && (mx||nx)) && (!(mx&&nx)))) {
295:       *inout = -*inout;
296:       continue;
297:     }
298:     val = (y[i]*x[j]-x[i]*y[j])/(x[j]-x[i]);
299:     if (PetscAbs(val)<10*PETSC_MACHINE_EPSILON) {
300:       *inout = 0;
301:       return(0);
302:     } else if (val>0.0) *inout = -*inout;
303:   }
304:   return(0);
305: }

309: PetscErrorCode RGSetFromOptions_Polygon(PetscOptionItems *PetscOptionsObject,RG rg)
310: {
312:   PetscScalar    array[VERTMAX];
313:   PetscInt       i,k;
314:   PetscBool      flg,flgi=PETSC_FALSE;
315: #if !defined(PETSC_USE_COMPLEX)
316:   PetscScalar    arrayi[VERTMAX];
317:   PetscInt       ki;
318: #else
319:   PetscScalar    *arrayi=NULL;
320: #endif

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

325:   k = VERTMAX;
326:   for (i=0;i<k;i++) array[i] = 0;
327:   PetscOptionsScalarArray("-rg_polygon_vertices","Vertices of polygon","RGPolygonSetVertices",array,&k,&flg);
328: #if !defined(PETSC_USE_COMPLEX)
329:   ki = VERTMAX;
330:   for (i=0;i<ki;i++) arrayi[i] = 0;
331:   PetscOptionsScalarArray("-rg_polygon_verticesi","Vertices of polygon (imaginary part)","RGPolygonSetVertices",arrayi,&ki,&flgi);
332:   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);
333: #endif
334:   if (flg || flgi) {
335:     RGPolygonSetVertices(rg,k,array,arrayi);
336:   }

338:   PetscOptionsTail();
339:   return(0);
340: }

344: PetscErrorCode RGDestroy_Polygon(RG rg)
345: {
347:   RG_POLYGON     *ctx = (RG_POLYGON*)rg->data;

350:   if (ctx->n) {
351:     PetscFree(ctx->vr);
352: #if !defined(PETSC_USE_COMPLEX)
353:     PetscFree(ctx->vi);
354: #endif
355:   }
356:   PetscFree(rg->data);
357:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",NULL);
358:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",NULL);
359:   return(0);
360: }

364: PETSC_EXTERN PetscErrorCode RGCreate_Polygon(RG rg)
365: {
366:   RG_POLYGON     *polygon;

370:   PetscNewLog(rg,&polygon);
371:   rg->data = (void*)polygon;

373:   rg->ops->istrivial      = RGIsTrivial_Polygon;
374:   rg->ops->computecontour = RGComputeContour_Polygon;
375:   rg->ops->checkinside    = RGCheckInside_Polygon;
376:   rg->ops->setfromoptions = RGSetFromOptions_Polygon;
377:   rg->ops->view           = RGView_Polygon;
378:   rg->ops->destroy        = RGDestroy_Polygon;
379:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonSetVertices_C",RGPolygonSetVertices_Polygon);
380:   PetscObjectComposeFunction((PetscObject)rg,"RGPolygonGetVertices_C",RGPolygonGetVertices_Polygon);
381:   return(0);
382: }