Actual source code: rgbasic.c

slepc-3.9.0 2018-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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:    Basic RG routines
 12: */

 14: #include <slepc/private/rgimpl.h>      /*I "slepcrg.h" I*/

 16: PetscFunctionList RGList = 0;
 17: PetscBool         RGRegisterAllCalled = PETSC_FALSE;
 18: PetscClassId      RG_CLASSID = 0;
 19: static PetscBool  RGPackageInitialized = PETSC_FALSE;

 21: /*@C
 22:    RGFinalizePackage - This function destroys everything in the Slepc interface
 23:    to the RG package. It is called from SlepcFinalize().

 25:    Level: developer

 27: .seealso: SlepcFinalize()
 28: @*/
 29: PetscErrorCode RGFinalizePackage(void)
 30: {

 34:   PetscFunctionListDestroy(&RGList);
 35:   RGPackageInitialized = PETSC_FALSE;
 36:   RGRegisterAllCalled  = PETSC_FALSE;
 37:   return(0);
 38: }

 40: /*@C
 41:   RGInitializePackage - This function initializes everything in the RG package.
 42:   It is called from PetscDLLibraryRegister() when using dynamic libraries, and
 43:   on the first call to RGCreate() when using static libraries.

 45:   Level: developer

 47: .seealso: SlepcInitialize()
 48: @*/
 49: PetscErrorCode RGInitializePackage(void)
 50: {
 51:   char           logList[256];
 52:   PetscBool      opt,pkg;

 56:   if (RGPackageInitialized) return(0);
 57:   RGPackageInitialized = PETSC_TRUE;
 58:   /* Register Classes */
 59:   PetscClassIdRegister("Region",&RG_CLASSID);
 60:   /* Register Constructors */
 61:   RGRegisterAll();
 62:   /* Process info exclusions */
 63:   PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,sizeof(logList),&opt);
 64:   if (opt) {
 65:     PetscStrInList("rg",logList,',',&pkg);
 66:     if (pkg) { PetscInfoDeactivateClass(RG_CLASSID); }
 67:   }
 68:   /* Process summary exclusions */
 69:   PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);
 70:   if (opt) {
 71:     PetscStrInList("rg",logList,',',&pkg);
 72:     if (pkg) { PetscLogEventDeactivateClass(RG_CLASSID); }
 73:   }
 74:   /* Register package finalizer */
 75:   PetscRegisterFinalize(RGFinalizePackage);
 76:   return(0);
 77: }

 79: /*@
 80:    RGCreate - Creates an RG context.

 82:    Collective on MPI_Comm

 84:    Input Parameter:
 85: .  comm - MPI communicator

 87:    Output Parameter:
 88: .  newrg - location to put the RG context

 90:    Level: beginner

 92: .seealso: RGDestroy(), RG
 93: @*/
 94: PetscErrorCode RGCreate(MPI_Comm comm,RG *newrg)
 95: {
 96:   RG             rg;

101:   *newrg = 0;
102:   RGInitializePackage();
103:   SlepcHeaderCreate(rg,RG_CLASSID,"RG","Region","RG",comm,RGDestroy,RGView);
104:   rg->complement = PETSC_FALSE;
105:   rg->sfactor    = 1.0;
106:   rg->osfactor   = 0.0;
107:   rg->data       = NULL;

109:   *newrg = rg;
110:   return(0);
111: }

113: /*@C
114:    RGSetOptionsPrefix - Sets the prefix used for searching for all
115:    RG options in the database.

117:    Logically Collective on RG

119:    Input Parameters:
120: +  rg     - the region context
121: -  prefix - the prefix string to prepend to all RG option requests

123:    Notes:
124:    A hyphen (-) must NOT be given at the beginning of the prefix name.
125:    The first character of all runtime options is AUTOMATICALLY the
126:    hyphen.

128:    Level: advanced

130: .seealso: RGAppendOptionsPrefix()
131: @*/
132: PetscErrorCode RGSetOptionsPrefix(RG rg,const char *prefix)
133: {

138:   PetscObjectSetOptionsPrefix((PetscObject)rg,prefix);
139:   return(0);
140: }

142: /*@C
143:    RGAppendOptionsPrefix - Appends to the prefix used for searching for all
144:    RG options in the database.

146:    Logically Collective on RG

148:    Input Parameters:
149: +  rg     - the region context
150: -  prefix - the prefix string to prepend to all RG option requests

152:    Notes:
153:    A hyphen (-) must NOT be given at the beginning of the prefix name.
154:    The first character of all runtime options is AUTOMATICALLY the hyphen.

156:    Level: advanced

158: .seealso: RGSetOptionsPrefix()
159: @*/
160: PetscErrorCode RGAppendOptionsPrefix(RG rg,const char *prefix)
161: {

166:   PetscObjectAppendOptionsPrefix((PetscObject)rg,prefix);
167:   return(0);
168: }

170: /*@C
171:    RGGetOptionsPrefix - Gets the prefix used for searching for all
172:    RG options in the database.

174:    Not Collective

176:    Input Parameters:
177: .  rg - the region context

179:    Output Parameters:
180: .  prefix - pointer to the prefix string used is returned

182:    Note:
183:    On the Fortran side, the user should pass in a string 'prefix' of
184:    sufficient length to hold the prefix.

186:    Level: advanced

188: .seealso: RGSetOptionsPrefix(), RGAppendOptionsPrefix()
189: @*/
190: PetscErrorCode RGGetOptionsPrefix(RG rg,const char *prefix[])
191: {

197:   PetscObjectGetOptionsPrefix((PetscObject)rg,prefix);
198:   return(0);
199: }

201: /*@C
202:    RGSetType - Selects the type for the RG object.

204:    Logically Collective on RG

206:    Input Parameter:
207: +  rg   - the region context
208: -  type - a known type

210:    Level: intermediate

212: .seealso: RGGetType()
213: @*/
214: PetscErrorCode RGSetType(RG rg,RGType type)
215: {
216:   PetscErrorCode ierr,(*r)(RG);
217:   PetscBool      match;


223:   PetscObjectTypeCompare((PetscObject)rg,type,&match);
224:   if (match) return(0);

226:    PetscFunctionListFind(RGList,type,&r);
227:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested RG type %s",type);

229:   if (rg->ops->destroy) { (*rg->ops->destroy)(rg); }
230:   PetscMemzero(rg->ops,sizeof(struct _RGOps));

232:   PetscObjectChangeTypeName((PetscObject)rg,type);
233:   (*r)(rg);
234:   return(0);
235: }

237: /*@C
238:    RGGetType - Gets the RG type name (as a string) from the RG context.

240:    Not Collective

242:    Input Parameter:
243: .  rg - the region context

245:    Output Parameter:
246: .  name - name of the region

248:    Level: intermediate

250: .seealso: RGSetType()
251: @*/
252: PetscErrorCode RGGetType(RG rg,RGType *type)
253: {
257:   *type = ((PetscObject)rg)->type_name;
258:   return(0);
259: }

261: /*@
262:    RGSetFromOptions - Sets RG options from the options database.

264:    Collective on RG

266:    Input Parameters:
267: .  rg - the region context

269:    Notes:
270:    To see all options, run your program with the -help option.

272:    Level: beginner
273: @*/
274: PetscErrorCode RGSetFromOptions(RG rg)
275: {
277:   char           type[256];
278:   PetscBool      flg;
279:   PetscReal      sfactor;

283:   RGRegisterAll();
284:   PetscObjectOptionsBegin((PetscObject)rg);
285:     PetscOptionsFList("-rg_type","Region type","RGSetType",RGList,(char*)(((PetscObject)rg)->type_name?((PetscObject)rg)->type_name:RGINTERVAL),type,256,&flg);
286:     if (flg) {
287:       RGSetType(rg,type);
288:     } else if (!((PetscObject)rg)->type_name) {
289:       RGSetType(rg,RGINTERVAL);
290:     }

292:     PetscOptionsBool("-rg_complement","Whether region is complemented or not","RGSetComplement",rg->complement,&rg->complement,NULL);

294:     PetscOptionsReal("-rg_scale","Scaling factor","RGSetScale",1.0,&sfactor,&flg);
295:     if (flg) { RGSetScale(rg,sfactor); }

297:     if (rg->ops->setfromoptions) {
298:       (*rg->ops->setfromoptions)(PetscOptionsObject,rg);
299:     }
300:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)rg);
301:   PetscOptionsEnd();
302:   return(0);
303: }

305: /*@C
306:    RGView - Prints the RG data structure.

308:    Collective on RG

310:    Input Parameters:
311: +  rg - the region context
312: -  viewer - optional visualization context

314:    Note:
315:    The available visualization contexts include
316: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
317: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
318:          output where only the first processor opens
319:          the file.  All other processors send their
320:          data to the first processor to print.

322:    The user can open an alternative visualization context with
323:    PetscViewerASCIIOpen() - output to a specified file.

325:    Level: beginner
326: @*/
327: PetscErrorCode RGView(RG rg,PetscViewer viewer)
328: {
329:   PetscBool      isdraw,isascii;
330:   PetscInt       tabs;

335:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)rg));
338:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
339:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
340:   if (isascii) {
341:     PetscViewerASCIIGetTab(viewer,&tabs);
342:     PetscViewerASCIISetTab(viewer,((PetscObject)rg)->tablevel);
343:     PetscObjectPrintClassNamePrefixType((PetscObject)rg,viewer);
344:     if (rg->ops->view) {
345:       PetscViewerASCIIPushTab(viewer);
346:       (*rg->ops->view)(rg,viewer);
347:       PetscViewerASCIIPopTab(viewer);
348:     }
349:     if (rg->complement) {
350:       PetscViewerASCIIPrintf(viewer,"  selected region is the complement of the specified one\n");
351:     }
352:     if (rg->sfactor!=1.0) {
353:       PetscViewerASCIIPrintf(viewer,"  scaling factor = %g\n",(double)rg->sfactor);
354:     }
355:     PetscViewerASCIISetTab(viewer,tabs);
356:   } else if (isdraw) {
357:     if (rg->ops->view) { (*rg->ops->view)(rg,viewer); }
358:   }
359:   return(0);
360: }

362: /*@
363:    RGIsTrivial - Whether it is the trivial region (whole complex plane).

365:    Not Collective

367:    Input Parameter:
368: .  rg - the region context

370:    Output Parameter:
371: .  trivial - true if the region is equal to the whole complex plane, e.g.,
372:              an interval region with all four endpoints unbounded or an
373:              ellipse with infinite radius.

375:    Level: beginner
376: @*/
377: PetscErrorCode RGIsTrivial(RG rg,PetscBool *trivial)
378: {

385:   if (rg->ops->istrivial) {
386:     (*rg->ops->istrivial)(rg,trivial);
387:   } else *trivial = PETSC_FALSE;
388:   return(0);
389: }

391: /*@
392:    RGCheckInside - Determines if a set of given points are inside the region or not.

394:    Not Collective

396:    Input Parameters:
397: +  rg - the region context
398: .  n  - number of points to check
399: .  ar - array of real parts
400: -  ai - array of imaginary parts

402:    Output Parameter:
403: .  inside - array of results (1=inside, 0=on the contour, -1=outside)

405:    Note:
406:    The point a is expressed as a couple of PetscScalar variables ar,ai.
407:    If built with complex scalars, the point is supposed to be stored in ar,
408:    otherwise ar,ai contain the real and imaginary parts, respectively.

410:    If a scaling factor was set, the points are scaled before checking.

412:    Level: intermediate

414: .seealso: RGSetScale(), RGSetComplement()
415: @*/
416: PetscErrorCode RGCheckInside(RG rg,PetscInt n,PetscScalar *ar,PetscScalar *ai,PetscInt *inside)
417: {
419:   PetscReal      px,py;
420:   PetscInt       i;

426: #if !defined(PETSC_USE_COMPLEX)
428: #endif

431:   for (i=0;i<n;i++) {
432: #if defined(PETSC_USE_COMPLEX)
433:     px = PetscRealPart(ar[i]);
434:     py = PetscImaginaryPart(ar[i]);
435: #else
436:     px = ar[i];
437:     py = ai[i];
438: #endif
439:     if (rg->sfactor != 1.0) {
440:       px /= rg->sfactor;
441:       py /= rg->sfactor;
442:     }
443:     (*rg->ops->checkinside)(rg,px,py,inside+i);
444:     if (rg->complement) inside[i] = -inside[i];
445:   }
446:   return(0);
447: }

449: /*@
450:    RGComputeContour - Computes the coordinates of several points lying in the
451:    contour of the region.

453:    Not Collective

455:    Input Parameters:
456: +  rg - the region context
457: -  n  - number of points to compute

459:    Output Parameters:
460: +  cr - location to store real parts
461: -  ci - location to store imaginary parts

463:    Level: intermediate
464: @*/
465: PetscErrorCode RGComputeContour(RG rg,PetscInt n,PetscScalar cr[],PetscScalar ci[])
466: {
467:   PetscInt       i;

474: #if !defined(PETSC_USE_COMPLEX)
476: #endif
477:   if (rg->complement) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Cannot compute contour of region with complement flag set");
478:   (*rg->ops->computecontour)(rg,n,cr,ci);
479:   for (i=0;i<n;i++) {
480:     cr[i] *= rg->sfactor;
481:     ci[i] *= rg->sfactor;
482:   }
483:   return(0);
484: }

486: /*@
487:    RGComputeBoundingBox - Determines the endpoints of a rectangle in the complex plane that
488:    contains the region.

490:    Not Collective

492:    Input Parameters:
493: .  rg - the region context

495:    Output Parameters:
496: +  a,b - endpoints of the bounding box in the real axis
497: -  c,d - endpoints of the bounding box in the imaginary axis

499:    Notes:
500:    The bounding box is defined as [a,b]x[c,d]. In regions that are not bounded (e.g. an
501:    open interval) or with the complement flag set, it makes no sense to compute a bounding
502:    box, so the return values are infinite.

504:    Level: intermediate

506: .seealso: RGSetComplement()
507: @*/
508: PetscErrorCode RGComputeBoundingBox(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
509: {


516:   if (rg->complement) {  /* cannot compute bounding box */
517:     if (a) *a = -PETSC_MAX_REAL;
518:     if (b) *b =  PETSC_MAX_REAL;
519:     if (c) *c = -PETSC_MAX_REAL;
520:     if (d) *d =  PETSC_MAX_REAL;
521:   } else {
522:     (*rg->ops->computebbox)(rg,a,b,c,d);
523:     if (a && *a!=-PETSC_MAX_REAL) *a *= rg->sfactor;
524:     if (b && *b!= PETSC_MAX_REAL) *b *= rg->sfactor;
525:     if (c && *c!=-PETSC_MAX_REAL) *c *= rg->sfactor;
526:     if (d && *d!= PETSC_MAX_REAL) *d *= rg->sfactor;
527:   }
528:   return(0);
529: }

531: /*@
532:    RGSetComplement - Sets a flag to indicate that the region is the complement
533:    of the specified one.

535:    Logically Collective on RG

537:    Input Parameters:
538: +  rg  - the region context
539: -  flg - the boolean flag

541:    Options Database Key:
542: .  -rg_complement <bool> - Activate/deactivate the complementation of the region

544:    Level: intermediate

546: .seealso: RGGetComplement()
547: @*/
548: PetscErrorCode RGSetComplement(RG rg,PetscBool flg)
549: {
553:   rg->complement = flg;
554:   return(0);
555: }

557: /*@
558:    RGGetComplement - Gets a flag that that indicates whether the region
559:    is complemented or not.

561:    Not Collective

563:    Input Parameter:
564: .  rg - the region context

566:    Output Parameter:
567: .  flg - the flag

569:    Level: intermediate

571: .seealso: RGSetComplement()
572: @*/
573: PetscErrorCode RGGetComplement(RG rg,PetscBool *flg)
574: {
578:   *flg = rg->complement;
579:   return(0);
580: }

582: /*@
583:    RGSetScale - Sets the scaling factor to be used when checking that a
584:    point is inside the region and when computing the contour.

586:    Logically Collective on RG

588:    Input Parameters:
589: +  rg      - the region context
590: -  sfactor - the scaling factor

592:    Options Database Key:
593: .  -rg_scale <real> - Sets the scaling factor

595:    Level: advanced

597: .seealso: RGGetScale(), RGCheckInside()
598: @*/
599: PetscErrorCode RGSetScale(RG rg,PetscReal sfactor)
600: {
604:   if (sfactor == PETSC_DEFAULT || sfactor == PETSC_DECIDE) rg->sfactor = 1.0;
605:   else {
606:     if (sfactor<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of scaling factor. Must be > 0");
607:     rg->sfactor = sfactor;
608:   }
609:   return(0);
610: }

612: /*@
613:    RGGetScale - Gets the scaling factor.

615:    Not Collective

617:    Input Parameter:
618: .  rg - the region context

620:    Output Parameter:
621: .  flg - the flag

623:    Level: advanced

625: .seealso: RGSetScale()
626: @*/
627: PetscErrorCode RGGetScale(RG rg,PetscReal *sfactor)
628: {
632:   *sfactor = rg->sfactor;
633:   return(0);
634: }

636: /*@
637:    RGPushScale - Sets an additional scaling factor, that will multiply the
638:    user-defined scaling factor.

640:    Logically Collective on RG

642:    Input Parameters:
643: +  rg      - the region context
644: -  sfactor - the scaling factor

646:    Notes:
647:    The current implementation does not allow pushing several scaling factors.

649:    This is intended for internal use, for instance in polynomial eigensolvers
650:    that use parameter scaling.

652:    Level: developer

654: .seealso: RGPopScale(), RGSetScale()
655: @*/
656: PetscErrorCode RGPushScale(RG rg,PetscReal sfactor)
657: {
661:   if (sfactor<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of scaling factor. Must be > 0");
662:   if (rg->osfactor) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Current implementation does not allow pushing several scaling factors");
663:   rg->osfactor = rg->sfactor;
664:   rg->sfactor *= sfactor;
665:   return(0);
666: }

668: /*@
669:    RGPopScale - Pops the scaling factor set with RGPushScale().

671:    Not Collective

673:    Input Parameter:
674: .  rg - the region context

676:    Level: developer

678: .seealso: RGPushScale()
679: @*/
680: PetscErrorCode RGPopScale(RG rg)
681: {
684:   if (!rg->osfactor) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ORDER,"Must call RGPushScale first");
685:   rg->sfactor  = rg->osfactor;
686:   rg->osfactor = 0.0;
687:   return(0);
688: }

690: /*@
691:    RGDestroy - Destroys RG context that was created with RGCreate().

693:    Collective on RG

695:    Input Parameter:
696: .  rg - the region context

698:    Level: beginner

700: .seealso: RGCreate()
701: @*/
702: PetscErrorCode RGDestroy(RG *rg)
703: {

707:   if (!*rg) return(0);
709:   if (--((PetscObject)(*rg))->refct > 0) { *rg = 0; return(0); }
710:   if ((*rg)->ops->destroy) { (*(*rg)->ops->destroy)(*rg); }
711:   PetscHeaderDestroy(rg);
712:   return(0);
713: }

715: /*@C
716:    RGRegister - Adds a region to the RG package.

718:    Not collective

720:    Input Parameters:
721: +  name - name of a new user-defined RG
722: -  function - routine to create context

724:    Notes:
725:    RGRegister() may be called multiple times to add several user-defined regions.

727:    Level: advanced

729: .seealso: RGRegisterAll()
730: @*/
731: PetscErrorCode RGRegister(const char *name,PetscErrorCode (*function)(RG))
732: {

736:   PetscFunctionListAdd(&RGList,name,function);
737:   return(0);
738: }