Actual source code: rginterval.c
slepc-3.7.1 2016-05-27
1: /*
2: Interval of the real axis or more generally a (possibly open) rectangle
3: of the complex plane.
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include <slepc/private/rgimpl.h> /*I "slepcrg.h" I*/
27: typedef struct {
28: PetscReal a,b; /* interval in the real axis */
29: PetscReal c,d; /* interval in the imaginary axis */
30: } RG_INTERVAL;
34: static PetscErrorCode RGIntervalSetEndpoints_Interval(RG rg,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
35: {
36: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
39: if (!a && !b && !c && !d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"At least one argument must be nonzero");
40: if (a==b && a) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, endpoints must be distinct (or both zero)");
41: if (a>b) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be a<b");
42: if (c==d && c) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, endpoints must be distinct (or both zero)");
43: if (c>d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be c<d");
44: #if !defined(PETSC_USE_COMPLEX)
45: if (c!=-d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"In real scalars the region must be symmetric wrt real axis");
46: #endif
47: ctx->a = a;
48: ctx->b = b;
49: ctx->c = c;
50: ctx->d = d;
51: return(0);
52: }
56: /*@
57: RGIntervalSetEndpoints - Sets the parameters defining the interval region.
59: Logically Collective on RG
61: Input Parameters:
62: + rg - the region context
63: . a,b - endpoints of the interval in the real axis
64: - c,d - endpoints of the interval in the imaginary axis
66: Options Database Keys:
67: . -rg_interval_endpoints - the four endpoints
69: Note:
70: The region is defined as [a,b]x[c,d]. Particular cases are an interval on
71: the real axis (c=d=0), similar for the imaginary axis (a=b=0), the whole
72: complex plane (a=-inf,b=inf,c=-inf,d=inf), and so on.
74: Level: advanced
76: .seealso: RGIntervalGetEndpoints()
77: @*/
78: PetscErrorCode RGIntervalSetEndpoints(RG rg,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
79: {
88: PetscTryMethod(rg,"RGIntervalSetEndpoints_C",(RG,PetscReal,PetscReal,PetscReal,PetscReal),(rg,a,b,c,d));
89: return(0);
90: }
94: static PetscErrorCode RGIntervalGetEndpoints_Interval(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
95: {
96: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
99: if (a) *a = ctx->a;
100: if (b) *b = ctx->b;
101: if (c) *c = ctx->c;
102: if (d) *d = ctx->d;
103: return(0);
104: }
108: /*@
109: RGIntervalGetEndpoints - Gets the parameters that define the interval region.
111: Not Collective
113: Input Parameter:
114: . rg - the region context
116: Output Parameters:
117: + a,b - endpoints of the interval in the real axis
118: - c,d - endpoints of the interval in the imaginary axis
120: Level: advanced
122: .seealso: RGIntervalSetEndpoints()
123: @*/
124: PetscErrorCode RGIntervalGetEndpoints(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
125: {
130: PetscUseMethod(rg,"RGIntervalGetEndpoints_C",(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(rg,a,b,c,d));
131: return(0);
132: }
136: PetscErrorCode RGView_Interval(RG rg,PetscViewer viewer)
137: {
139: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
140: PetscBool isascii;
143: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
144: if (isascii) {
145: PetscViewerASCIIPrintf(viewer,"region: [%g,%g]x[%g,%g]\n",RGShowReal(ctx->a),RGShowReal(ctx->b),RGShowReal(ctx->c),RGShowReal(ctx->d));
146: }
147: return(0);
148: }
152: PetscErrorCode RGIsTrivial_Interval(RG rg,PetscBool *trivial)
153: {
154: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
157: if (rg->complement) *trivial = (ctx->a==ctx->b && ctx->c==ctx->d)? PETSC_TRUE: PETSC_FALSE;
158: 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;
159: return(0);
160: }
164: PetscErrorCode RGComputeContour_Interval(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
165: {
166: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
167: PetscInt i,pt,idx,j;
168: PetscReal hr[4],hi[4],h,off,d[4];
169: PetscScalar vr[4],vi[4];
172: 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");
173: if (ctx->a==ctx->b || ctx->c==ctx->d) {
174: if (ctx->a==ctx->b) {hi[0] = (ctx->d-ctx->c)/(n-1); hr[0] = 0.0;}
175: else {hr[0] = (ctx->b-ctx->a)/(n-1); hi[0] = 0.0;}
176: for (i=0;i<n;i++) {
177: #if defined(PETSC_USE_COMPLEX)
178: cr[i] = ctx->a+hr[0]*i + (ctx->c+hi[0]*i)*PETSC_i;
179: #else
180: cr[i] = ctx->a+hr[0]*i; ci[i] = ctx->c+hi[0]*i;
181: #endif
182: }
183: } else {
184: d[1] = d[3] = ctx->d-ctx->c; d[0] = d[2] = ctx->b-ctx->a;
185: h = (2*(d[0]+d[1]))/n;
186: vr[0] = ctx->a; vr[1] = ctx->b; vr[2] = ctx->b; vr[3] = ctx->a;
187: vi[0] = ctx->c; vi[1] = ctx->c; vi[2] = ctx->d; vi[3] = ctx->d;
188: hr[0] = h; hr[1] = 0.0; hr[2] = -h; hr[3] = 0.0;
189: hi[0] = 0.0; hi[1] = h; hi[2] = 0.0; hi[3] = -h;
190: off = 0.0; idx = 0;
191: for (i=0;i<4;i++) {
192: pt = (d[i]-off)/h+1;
193: #if defined(PETSC_USE_COMPLEX)
194: cr[idx] = vr[i]+off*(hr[i]/h)+ (vi[i]+off*(hi[i]/h))*PETSC_i;
195: #else
196: cr[idx] = vr[i]+off*(hr[i]/h); ci[idx]=vi[i]+off*(hi[i]/h);
197: #endif
198: idx++;
199: for (j=1;j<pt;j++) {
200: #if defined(PETSC_USE_COMPLEX)
201: cr[idx] = cr[idx-1]+(hr[i]+hi[i]*PETSC_i);
202: #else
203: cr[idx] = cr[idx-1]+hr[i]; ci[idx] = ci[idx-1]+hi[i];
204: #endif
205: idx++;
206: }
207: off = off+pt*h-d[i];
208: if (off>=d[i+1]) {off -= d[i+1]; i++;}
209: }
210: }
211: return(0);
212: }
216: PetscErrorCode RGCheckInside_Interval(RG rg,PetscReal dx,PetscReal dy,PetscInt *inside)
217: {
218: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
221: if (dx>ctx->a && dx<ctx->b) *inside = 1;
222: else if (dx==ctx->a || dx==ctx->b) *inside = 0;
223: else *inside = -1;
224: if (*inside>=0) {
225: if (dy>ctx->c && dy<ctx->d) ;
226: else if (dy==ctx->c || dy==ctx->d) *inside = 0;
227: else *inside = -1;
228: }
229: return(0);
230: }
234: PetscErrorCode RGSetFromOptions_Interval(PetscOptionItems *PetscOptionsObject,RG rg)
235: {
237: PetscBool flg;
238: PetscInt k;
239: PetscReal array[4]={0,0,0,0};
242: PetscOptionsHead(PetscOptionsObject,"RG Interval Options");
244: k = 4;
245: PetscOptionsRealArray("-rg_interval_endpoints","Interval endpoints (four real values separated with a comma without spaces)","RGIntervalSetEndpoints",array,&k,&flg);
246: if (flg) {
247: if (k<2) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_SIZ,"Must pass at leat two values in -rg_interval_endpoints (comma-separated without spaces)");
248: RGIntervalSetEndpoints(rg,array[0],array[1],array[2],array[3]);
249: }
251: PetscOptionsTail();
252: return(0);
253: }
257: PetscErrorCode RGDestroy_Interval(RG rg)
258: {
262: PetscFree(rg->data);
263: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",NULL);
264: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",NULL);
265: return(0);
266: }
270: PETSC_EXTERN PetscErrorCode RGCreate_Interval(RG rg)
271: {
272: RG_INTERVAL *interval;
276: PetscNewLog(rg,&interval);
277: interval->a = -PETSC_MAX_REAL;
278: interval->b = PETSC_MAX_REAL;
279: interval->c = -PETSC_MAX_REAL;
280: interval->d = PETSC_MAX_REAL;
281: rg->data = (void*)interval;
283: rg->ops->istrivial = RGIsTrivial_Interval;
284: rg->ops->computecontour = RGComputeContour_Interval;
285: rg->ops->checkinside = RGCheckInside_Interval;
286: rg->ops->setfromoptions = RGSetFromOptions_Interval;
287: rg->ops->view = RGView_Interval;
288: rg->ops->destroy = RGDestroy_Interval;
289: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",RGIntervalSetEndpoints_Interval);
290: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",RGIntervalGetEndpoints_Interval);
291: return(0);
292: }