Actual source code: rginterval.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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-2014, 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:   PetscTryMethod(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;

169:   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");
170:   SETERRQ(PetscObjectComm((PetscObject)rg),1,"Not implemented yet");
171:   return(0);
172: }

176: PetscErrorCode RGCheckInside_Interval(RG rg,PetscInt n,PetscScalar *ar,PetscScalar *ai,PetscInt *inside)
177: {
178:   RG_INTERVAL *ctx = (RG_INTERVAL*)rg->data;
179:   PetscInt    i;
180:   PetscReal   dx,dy;

183:   for (i=0;i<n;i++) {
184: #if defined(PETSC_USE_COMPLEX)
185:     dx = PetscRealPart(ar[i]);
186:     dy = PetscImaginaryPart(ar[i]);
187: #else
188:     dx = ar[i];
189:     dy = ai[i];
190: #endif
191:     if (dx>ctx->a && dx<ctx->b) inside[i] = 1;
192:     else if (dx==ctx->a || dx==ctx->b) inside[i] = 0;
193:     else inside[i] = -1;
194:     if (inside[i]>=0) {
195:       if (dy>ctx->c && dy<ctx->d) ;
196:       else if (dy==ctx->c || dy==ctx->d) inside[i] = 0;
197:       else inside[i] = -1;
198:     }
199:   }
200:   return(0);
201: }

205: PetscErrorCode RGSetFromOptions_Interval(RG rg)
206: {
208:   PetscBool      flg;
209:   PetscInt       k;
210:   PetscReal      array[4]={0,0,0,0};

213:   PetscOptionsHead("RG Interval Options");

215:   k = 4;
216:   PetscOptionsRealArray("-rg_interval_endpoints","Interval endpoints (four real values separated with a comma without spaces)","RGIntervalSetEndpoints",array,&k,&flg);
217:   if (flg) {
218:     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)");
219:     RGIntervalSetEndpoints(rg,array[0],array[1],array[2],array[3]);
220:   }

222:   PetscOptionsTail();
223:   return(0);
224: }

228: PetscErrorCode RGDestroy_Interval(RG rg)
229: {

233:   PetscFree(rg->data);
234:   PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",NULL);
235:   PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",NULL);
236:   return(0);
237: }

241: PETSC_EXTERN PetscErrorCode RGCreate_Interval(RG rg)
242: {
243:   RG_INTERVAL    *interval;

247:   PetscNewLog(rg,&interval);
248:   interval->a = -PETSC_MAX_REAL;
249:   interval->b = PETSC_MAX_REAL;
250:   interval->c = -PETSC_MAX_REAL;
251:   interval->d = PETSC_MAX_REAL;
252:   rg->data = (void*)interval;

254:   rg->ops->istrivial      = RGIsTrivial_Interval;
255:   rg->ops->computecontour = RGComputeContour_Interval;
256:   rg->ops->checkinside    = RGCheckInside_Interval;
257:   rg->ops->setfromoptions = RGSetFromOptions_Interval;
258:   rg->ops->view           = RGView_Interval;
259:   rg->ops->destroy        = RGDestroy_Interval;
260:   PetscObjectComposeFunction((PetscObject)rg,"RGIntervalSetEndpoints_C",RGIntervalSetEndpoints_Interval);
261:   PetscObjectComposeFunction((PetscObject)rg,"RGIntervalGetEndpoints_C",RGIntervalGetEndpoints_Interval);
262:   return(0);
263: }