Actual source code: rginterval.c
slepc-3.8.2 2017-12-01
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: Interval of the real axis or more generally a (possibly open) rectangle
12: of the complex plane
13: */
15: #include <slepc/private/rgimpl.h> /*I "slepcrg.h" I*/
17: typedef struct {
18: PetscReal a,b; /* interval in the real axis */
19: PetscReal c,d; /* interval in the imaginary axis */
20: } RG_INTERVAL;
22: static PetscErrorCode RGIntervalSetEndpoints_Interval(RG rg,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
23: {
24: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
27: if (!a && !b && !c && !d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"At least one argument must be nonzero");
28: if (a==b && a) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, endpoints must be distinct (or both zero)");
29: if (a>b) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be a<b");
30: if (c==d && c) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, endpoints must be distinct (or both zero)");
31: if (c>d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be c<d");
32: #if !defined(PETSC_USE_COMPLEX)
33: if (c!=-d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"In real scalars the region must be symmetric wrt real axis");
34: #endif
35: ctx->a = a;
36: ctx->b = b;
37: ctx->c = c;
38: ctx->d = d;
39: return(0);
40: }
42: /*@
43: RGIntervalSetEndpoints - Sets the parameters defining the interval region.
45: Logically Collective on RG
47: Input Parameters:
48: + rg - the region context
49: . a,b - endpoints of the interval in the real axis
50: - c,d - endpoints of the interval in the imaginary axis
52: Options Database Keys:
53: . -rg_interval_endpoints - the four endpoints
55: Note:
56: The region is defined as [a,b]x[c,d]. Particular cases are an interval on
57: the real axis (c=d=0), similar for the imaginary axis (a=b=0), the whole
58: complex plane (a=-inf,b=inf,c=-inf,d=inf), and so on.
60: When PETSc is built with real scalars, the region must be symmetric with
61: respect to the real axis.
63: Level: advanced
65: .seealso: RGIntervalGetEndpoints()
66: @*/
67: PetscErrorCode RGIntervalSetEndpoints(RG rg,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
68: {
77: PetscTryMethod(rg,"RGIntervalSetEndpoints_C",(RG,PetscReal,PetscReal,PetscReal,PetscReal),(rg,a,b,c,d));
78: return(0);
79: }
81: static PetscErrorCode RGIntervalGetEndpoints_Interval(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
82: {
83: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
86: if (a) *a = ctx->a;
87: if (b) *b = ctx->b;
88: if (c) *c = ctx->c;
89: if (d) *d = ctx->d;
90: return(0);
91: }
93: /*@
94: RGIntervalGetEndpoints - Gets the parameters that define the interval region.
96: Not Collective
98: Input Parameter:
99: . rg - the region context
101: Output Parameters:
102: + a,b - endpoints of the interval in the real axis
103: - c,d - endpoints of the interval in the imaginary axis
105: Level: advanced
107: .seealso: RGIntervalSetEndpoints()
108: @*/
109: PetscErrorCode RGIntervalGetEndpoints(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
110: {
115: PetscUseMethod(rg,"RGIntervalGetEndpoints_C",(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(rg,a,b,c,d));
116: return(0);
117: }
119: PetscErrorCode RGView_Interval(RG rg,PetscViewer viewer)
120: {
122: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
123: PetscBool isascii;
126: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
127: if (isascii) {
128: PetscViewerASCIIPrintf(viewer," region: [%g,%g]x[%g,%g]\n",RGShowReal(ctx->a),RGShowReal(ctx->b),RGShowReal(ctx->c),RGShowReal(ctx->d));
129: }
130: return(0);
131: }
133: PetscErrorCode RGIsTrivial_Interval(RG rg,PetscBool *trivial)
134: {
135: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
138: if (rg->complement) *trivial = (ctx->a==ctx->b && ctx->c==ctx->d)? PETSC_TRUE: PETSC_FALSE;
139: else *trivial = (ctx->a<=-PETSC_MAX_REAL && ctx->b>=PETSC_MAX_REAL && ctx->c<=-PETSC_MAX_REAL && ctx->d>=PETSC_MAX_REAL)? PETSC_TRUE: PETSC_FALSE;
140: return(0);
141: }
143: PetscErrorCode RGComputeContour_Interval(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
144: {
145: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
146: PetscInt i,pt,idx,j;
147: PetscReal hr[4],hi[4],h,off,d[4];
148: PetscScalar vr[4],vi[4];
151: if (!(ctx->a>-PETSC_MAX_REAL && ctx->b<PETSC_MAX_REAL && ctx->c>-PETSC_MAX_REAL && ctx->d<PETSC_MAX_REAL)) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Contour not defined in unbounded regions");
152: if (ctx->a==ctx->b || ctx->c==ctx->d) {
153: if (ctx->a==ctx->b) {hi[0] = (ctx->d-ctx->c)/(n-1); hr[0] = 0.0;}
154: else {hr[0] = (ctx->b-ctx->a)/(n-1); hi[0] = 0.0;}
155: for (i=0;i<n;i++) {
156: #if defined(PETSC_USE_COMPLEX)
157: cr[i] = ctx->a+hr[0]*i + (ctx->c+hi[0]*i)*PETSC_i;
158: #else
159: cr[i] = ctx->a+hr[0]*i; ci[i] = ctx->c+hi[0]*i;
160: #endif
161: }
162: } else {
163: d[1] = d[3] = ctx->d-ctx->c; d[0] = d[2] = ctx->b-ctx->a;
164: h = (2*(d[0]+d[1]))/n;
165: vr[0] = ctx->a; vr[1] = ctx->b; vr[2] = ctx->b; vr[3] = ctx->a;
166: vi[0] = ctx->c; vi[1] = ctx->c; vi[2] = ctx->d; vi[3] = ctx->d;
167: hr[0] = h; hr[1] = 0.0; hr[2] = -h; hr[3] = 0.0;
168: hi[0] = 0.0; hi[1] = h; hi[2] = 0.0; hi[3] = -h;
169: off = 0.0; idx = 0;
170: for (i=0;i<4;i++) {
171: #if defined(PETSC_USE_COMPLEX)
172: cr[idx] = vr[i]+off*(hr[i]/h)+ (vi[i]+off*(hi[i]/h))*PETSC_i;
173: #else
174: cr[idx] = vr[i]+off*(hr[i]/h); ci[idx]=vi[i]+off*(hi[i]/h);
175: #endif
176: idx++;
177: pt = (d[i]-off)/h+1;
178: for (j=1;j<pt && idx<n;j++) {
179: #if defined(PETSC_USE_COMPLEX)
180: cr[idx] = cr[idx-1]+(hr[i]+hi[i]*PETSC_i);
181: #else
182: cr[idx] = cr[idx-1]+hr[i]; ci[idx] = ci[idx-1]+hi[i];
183: #endif
184: idx++;
185: }
186: off = off+pt*h-d[i];
187: if (off>=d[i+1]) {off -= d[i+1]; i++;}
188: }
189: }
190: return(0);
191: }
193: PetscErrorCode RGComputeBoundingBox_Interval(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
194: {
195: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
198: if (a) *a = ctx->a;
199: if (b) *b = ctx->b;
200: if (c) *c = ctx->c;
201: if (d) *d = ctx->d;
202: return(0);
203: }
205: PetscErrorCode RGCheckInside_Interval(RG rg,PetscReal dx,PetscReal dy,PetscInt *inside)
206: {
207: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
210: if (dx>ctx->a && dx<ctx->b) *inside = 1;
211: else if (dx==ctx->a || dx==ctx->b) *inside = 0;
212: else *inside = -1;
213: if (*inside>=0) {
214: if (dy>ctx->c && dy<ctx->d) ;
215: else if (dy==ctx->c || dy==ctx->d) *inside = 0;
216: else *inside = -1;
217: }
218: return(0);
219: }
221: PetscErrorCode RGSetFromOptions_Interval(PetscOptionItems *PetscOptionsObject,RG rg)
222: {
224: PetscBool flg;
225: PetscInt k;
226: PetscReal array[4]={0,0,0,0};
229: PetscOptionsHead(PetscOptionsObject,"RG Interval Options");
231: k = 4;
232: PetscOptionsRealArray("-rg_interval_endpoints","Interval endpoints (two or four real values separated with a comma without spaces)","RGIntervalSetEndpoints",array,&k,&flg);
233: if (flg) {
234: if (k<2) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_SIZ,"Must pass at least two values in -rg_interval_endpoints (comma-separated without spaces)");
235: RGIntervalSetEndpoints(rg,array[0],array[1],array[2],array[3]);
236: }
238: PetscOptionsTail();
239: return(0);
240: }
242: PetscErrorCode RGDestroy_Interval(RG rg)
243: {
247: PetscFree(rg->data);
248: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",NULL);
249: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",NULL);
250: return(0);
251: }
253: PETSC_EXTERN PetscErrorCode RGCreate_Interval(RG rg)
254: {
255: RG_INTERVAL *interval;
259: PetscNewLog(rg,&interval);
260: interval->a = -PETSC_MAX_REAL;
261: interval->b = PETSC_MAX_REAL;
262: interval->c = -PETSC_MAX_REAL;
263: interval->d = PETSC_MAX_REAL;
264: rg->data = (void*)interval;
266: rg->ops->istrivial = RGIsTrivial_Interval;
267: rg->ops->computecontour = RGComputeContour_Interval;
268: rg->ops->computebbox = RGComputeBoundingBox_Interval;
269: rg->ops->checkinside = RGCheckInside_Interval;
270: rg->ops->setfromoptions = RGSetFromOptions_Interval;
271: rg->ops->view = RGView_Interval;
272: rg->ops->destroy = RGDestroy_Interval;
273: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",RGIntervalSetEndpoints_Interval);
274: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",RGIntervalGetEndpoints_Interval);
275: return(0);
276: }