Actual source code: rginterval.c
slepc-3.16.1 2021-11-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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>
16: #include <petscdraw.h>
18: typedef struct {
19: PetscReal a,b; /* interval in the real axis */
20: PetscReal c,d; /* interval in the imaginary axis */
21: } RG_INTERVAL;
23: static PetscErrorCode RGIntervalSetEndpoints_Interval(RG rg,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
24: {
25: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
28: if (!a && !b && !c && !d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"At least one argument must be nonzero");
29: if (a==b && a) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, endpoints must be distinct (or both zero)");
30: if (a>b) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be a<b");
31: if (c==d && c) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, endpoints must be distinct (or both zero)");
32: if (c>d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be c<d");
33: #if !defined(PETSC_USE_COMPLEX)
34: if (c!=-d) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"In real scalars the region must be symmetric wrt real axis");
35: #endif
36: ctx->a = a;
37: ctx->b = b;
38: ctx->c = c;
39: ctx->d = d;
40: return(0);
41: }
43: /*@
44: RGIntervalSetEndpoints - Sets the parameters defining the interval region.
46: Logically Collective on rg
48: Input Parameters:
49: + rg - the region context
50: . a - left endpoint of the interval in the real axis
51: . b - right endpoint of the interval in the real axis
52: . c - bottom endpoint of the interval in the imaginary axis
53: - d - top endpoint of the interval in the imaginary axis
55: Options Database Keys:
56: . -rg_interval_endpoints - the four endpoints
58: Note:
59: The region is defined as [a,b]x[c,d]. Particular cases are an interval on
60: the real axis (c=d=0), similar for the imaginary axis (a=b=0), the whole
61: complex plane (a=-inf,b=inf,c=-inf,d=inf), and so on.
63: When PETSc is built with real scalars, the region must be symmetric with
64: respect to the real axis.
66: Level: advanced
68: .seealso: RGIntervalGetEndpoints()
69: @*/
70: PetscErrorCode RGIntervalSetEndpoints(RG rg,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
71: {
80: PetscTryMethod(rg,"RGIntervalSetEndpoints_C",(RG,PetscReal,PetscReal,PetscReal,PetscReal),(rg,a,b,c,d));
81: return(0);
82: }
84: static PetscErrorCode RGIntervalGetEndpoints_Interval(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
85: {
86: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
89: if (a) *a = ctx->a;
90: if (b) *b = ctx->b;
91: if (c) *c = ctx->c;
92: if (d) *d = ctx->d;
93: return(0);
94: }
96: /*@C
97: RGIntervalGetEndpoints - Gets the parameters that define the interval region.
99: Not Collective
101: Input Parameter:
102: . rg - the region context
104: Output Parameters:
105: + a - left endpoint of the interval in the real axis
106: . b - right endpoint of the interval in the real axis
107: . c - bottom endpoint of the interval in the imaginary axis
108: - d - top endpoint of the interval in the imaginary axis
110: Level: advanced
112: .seealso: RGIntervalSetEndpoints()
113: @*/
114: PetscErrorCode RGIntervalGetEndpoints(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
115: {
120: PetscUseMethod(rg,"RGIntervalGetEndpoints_C",(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(rg,a,b,c,d));
121: return(0);
122: }
124: PetscErrorCode RGView_Interval(RG rg,PetscViewer viewer)
125: {
127: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
128: PetscBool isdraw,isascii;
129: int winw,winh;
130: PetscDraw draw;
131: PetscDrawAxis axis;
132: PetscReal a,b,c,d,ab,cd,lx,ly,w,scale=1.2;
135: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
136: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
137: if (isascii) {
138: PetscViewerASCIIPrintf(viewer," region: [%g,%g]x[%g,%g]\n",RGShowReal(ctx->a),RGShowReal(ctx->b),RGShowReal(ctx->c),RGShowReal(ctx->d));
139: } else if (isdraw) {
140: PetscViewerDrawGetDraw(viewer,0,&draw);
141: PetscDrawCheckResizedWindow(draw);
142: PetscDrawGetWindowSize(draw,&winw,&winh);
143: winw = PetscMax(winw,1); winh = PetscMax(winh,1);
144: PetscDrawClear(draw);
145: PetscDrawSetTitle(draw,"Interval region");
146: PetscDrawAxisCreate(draw,&axis);
147: a = ctx->a*rg->sfactor;
148: b = ctx->b*rg->sfactor;
149: c = ctx->c*rg->sfactor;
150: d = ctx->d*rg->sfactor;
151: lx = b-a;
152: ly = d-c;
153: ab = (a+b)/2;
154: cd = (c+d)/2;
155: w = scale*PetscMax(lx/winw,ly/winh)/2;
156: PetscDrawAxisSetLimits(axis,ab-w*winw,ab+w*winw,cd-w*winh,cd+w*winh);
157: PetscDrawAxisDraw(axis);
158: PetscDrawAxisDestroy(&axis);
159: PetscDrawRectangle(draw,a,c,b,d,PETSC_DRAW_BLUE,PETSC_DRAW_BLUE,PETSC_DRAW_BLUE,PETSC_DRAW_BLUE);
160: PetscDrawFlush(draw);
161: PetscDrawSave(draw);
162: PetscDrawPause(draw);
163: }
164: return(0);
165: }
167: PetscErrorCode RGIsTrivial_Interval(RG rg,PetscBool *trivial)
168: {
169: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
172: if (rg->complement) *trivial = (ctx->a==ctx->b && ctx->c==ctx->d)? PETSC_TRUE: PETSC_FALSE;
173: 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;
174: return(0);
175: }
177: PetscErrorCode RGComputeContour_Interval(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
178: {
179: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
180: PetscInt i,N,Nv,Nh,k1,k0;
181: PetscReal hv,hh,t;
184: 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");
185: if (ctx->a==ctx->b || ctx->c==ctx->d) {
186: if (ctx->a==ctx->b) {hv = (ctx->d-ctx->c)/(n-1); hh = 0.0;}
187: else {hh = (ctx->b-ctx->a)/(n-1); hv = 0.0;}
188: for (i=0;i<n;i++) {
189: #if defined(PETSC_USE_COMPLEX)
190: cr[i] = PetscCMPLX(ctx->a+hh*i,ctx->c+hv*i);
191: #else
192: if (cr) cr[i] = ctx->a+hh*i;
193: if (ci) ci[i] = ctx->c+hv*i;
194: #endif
195: }
196: } else {
197: if (n<4) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Minimum number of contour points: 4");
198: N = n/2;
199: t = ((ctx->d-ctx->c)/(ctx->d-ctx->c+ctx->b-ctx->a))*N;
200: Nv = t-PetscFloorReal(t)>0.5?PetscCeilReal(t):PetscFloorReal(t);
201: if (Nv==0) Nv++;
202: else if (Nv==N) Nv--;
203: Nh = N-Nv;
204: hh = (ctx->b-ctx->a)/Nh;
205: hv = (ctx->d-ctx->c)/Nv;
206: /* positive imaginary part first */
207: k1 = Nv/2+1;
208: k0 = Nv-k1;
210: for (i=k1;i<Nv;i++) {
211: #if defined(PETSC_USE_COMPLEX)
212: cr[i-k1] = PetscCMPLX(ctx->b,ctx->c+i*hv);
213: cr[i-k1+N] = PetscCMPLX(ctx->a,ctx->d-i*hv);
214: #else
215: if (cr) {cr[i-k1] = ctx->b; cr[i-k1+N] = ctx->a;}
216: if (ci) {ci[i-k1] = ctx->c+i*hv; ci[i-k1+N] = ctx->d-i*hv;}
217: #endif
218: }
219: for (i=0;i<Nh;i++) {
220: #if defined(PETSC_USE_COMPLEX)
221: cr[i+k0] = PetscCMPLX(ctx->b-i*hh,ctx->d);
222: cr[i+k0+N] = PetscCMPLX(ctx->a+i*hh,ctx->c);
223: #else
224: if (cr) {cr[i+k0] = ctx->b-i*hh; cr[i+k0+N] = ctx->a+i*hh;}
225: if (ci) {ci[i+k0] = ctx->d; ci[i+k0+N] = ctx->c;}
226: #endif
227: }
228: for (i=0;i<k1;i++) {
229: #if defined(PETSC_USE_COMPLEX)
230: cr[i+k0+Nh] = PetscCMPLX(ctx->a,ctx->d-i*hv);
231: cr[i+k0+Nh+N] = PetscCMPLX(ctx->b,ctx->c+i*hv);
232: #else
233: if (cr) {cr[i+k0+Nh] = ctx->a; cr[i+k0+Nh+N] = ctx->b;}
234: if (ci) {ci[i+k0+Nh] = ctx->d+i*hv; ci[i+k0+Nh+N] = ctx->c-i*hv;}
235: #endif
236: }
237: }
238: return(0);
239: }
241: PetscErrorCode RGComputeBoundingBox_Interval(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
242: {
243: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
246: if (a) *a = ctx->a;
247: if (b) *b = ctx->b;
248: if (c) *c = ctx->c;
249: if (d) *d = ctx->d;
250: return(0);
251: }
253: PetscErrorCode RGComputeQuadrature_Interval(RG rg,RGQuadRule quad,PetscInt n,PetscScalar *z,PetscScalar *zn,PetscScalar *w)
254: {
255: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
256: PetscReal theta,max_w=0.0,radius=1.0;
257: PetscScalar tmp,tmp2,center=0.0;
258: PetscInt i,j;
261: if (quad == RG_QUADRULE_CHEBYSHEV) {
262: if ((ctx->c!=ctx->d || ctx->c!=0.0) && (ctx->a!=ctx->b || ctx->a!=0.0)) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Endpoints of the imaginary axis or the real axis must be both zero");
263: for (i=0;i<n;i++) {
264: theta = PETSC_PI*(i+0.5)/n;
265: zn[i] = PetscCosReal(theta);
266: w[i] = PetscCosReal((n-1)*theta)/n;
267: if (ctx->c==ctx->d) z[i] = ((ctx->b-ctx->a)*(zn[i]+1.0)/2.0+ctx->a)*rg->sfactor;
268: else if (ctx->a==ctx->b) {
269: #if defined(PETSC_USE_COMPLEX)
270: z[i] = ((ctx->d-ctx->c)*(zn[i]+1.0)/2.0+ctx->c)*rg->sfactor*PETSC_i;
271: #else
272: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
273: #endif
274: }
275: }
276: } else { /* RG_QUADRULE_TRAPEZOIDAL */
277: #if defined(PETSC_USE_COMPLEX)
278: center = rg->sfactor*PetscCMPLX(ctx->b+ctx->a,ctx->d+ctx->c)/2.0;
279: #else
280: center = rg->sfactor*(ctx->b+ctx->a)/2.0;
281: #endif
282: radius = PetscSqrtReal(PetscPowRealInt(rg->sfactor*(ctx->b-ctx->a)/2.0,2)+PetscPowRealInt(rg->sfactor*(ctx->d-ctx->c)/2.0,2));
283: for (i=0;i<n;i++) {
284: zn[i] = (z[i]-center)/radius;
285: tmp = 1.0; tmp2 = 1.0;
286: for (j=0;j<n;j++) {
287: tmp *= z[j];
288: if (i != j) tmp2 *= z[j]-z[i];
289: }
290: w[i] = tmp/tmp2;
291: max_w = PetscMax(PetscAbsScalar(w[i]),max_w);
292: }
293: for (i=0;i<n;i++) w[i] /= (PetscScalar)max_w;
294: }
295: return(0);
296: }
298: PetscErrorCode RGCheckInside_Interval(RG rg,PetscReal dx,PetscReal dy,PetscInt *inside)
299: {
300: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
303: if (dx>ctx->a && dx<ctx->b) *inside = 1;
304: else if (dx==ctx->a || dx==ctx->b) *inside = 0;
305: else *inside = -1;
306: if (*inside>=0) {
307: if (dy>ctx->c && dy<ctx->d) ;
308: else if (dy==ctx->c || dy==ctx->d) *inside = 0;
309: else *inside = -1;
310: }
311: return(0);
312: }
314: PetscErrorCode RGIsAxisymmetric_Interval(RG rg,PetscBool vertical,PetscBool *symm)
315: {
316: RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
319: if (vertical) *symm = (ctx->a == -ctx->b)? PETSC_TRUE: PETSC_FALSE;
320: else *symm = (ctx->c == -ctx->d)? PETSC_TRUE: PETSC_FALSE;
321: return(0);
322: }
324: PetscErrorCode RGSetFromOptions_Interval(PetscOptionItems *PetscOptionsObject,RG rg)
325: {
327: PetscBool flg;
328: PetscInt k;
329: PetscReal array[4]={0,0,0,0};
332: PetscOptionsHead(PetscOptionsObject,"RG Interval Options");
334: k = 4;
335: PetscOptionsRealArray("-rg_interval_endpoints","Interval endpoints (two or four real values separated with a comma without spaces)","RGIntervalSetEndpoints",array,&k,&flg);
336: if (flg) {
337: 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)");
338: RGIntervalSetEndpoints(rg,array[0],array[1],array[2],array[3]);
339: }
341: PetscOptionsTail();
342: return(0);
343: }
345: PetscErrorCode RGDestroy_Interval(RG rg)
346: {
350: PetscFree(rg->data);
351: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",NULL);
352: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",NULL);
353: return(0);
354: }
356: SLEPC_EXTERN PetscErrorCode RGCreate_Interval(RG rg)
357: {
358: RG_INTERVAL *interval;
362: PetscNewLog(rg,&interval);
363: interval->a = -PETSC_MAX_REAL;
364: interval->b = PETSC_MAX_REAL;
365: interval->c = -PETSC_MAX_REAL;
366: interval->d = PETSC_MAX_REAL;
367: rg->data = (void*)interval;
369: rg->ops->istrivial = RGIsTrivial_Interval;
370: rg->ops->computecontour = RGComputeContour_Interval;
371: rg->ops->computebbox = RGComputeBoundingBox_Interval;
372: rg->ops->computequadrature = RGComputeQuadrature_Interval;
373: rg->ops->checkinside = RGCheckInside_Interval;
374: rg->ops->isaxisymmetric = RGIsAxisymmetric_Interval;
375: rg->ops->setfromoptions = RGSetFromOptions_Interval;
376: rg->ops->view = RGView_Interval;
377: rg->ops->destroy = RGDestroy_Interval;
378: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",RGIntervalSetEndpoints_Interval);
379: PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",RGIntervalGetEndpoints_Interval);
380: return(0);
381: }