Actual source code: rgbasic.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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;

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

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

363:    Not Collective

365:    Input Parameter:
366: .  rg - the region context

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

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

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

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

392:    Not Collective

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

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

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

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

410:    Level: intermediate

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

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

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

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

451:    Not Collective

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

457:    Output Parameters:
458: +  cr - location to store real parts
459: -  ci - location to store imaginary parts

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

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

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

488:    Not Collective

490:    Input Parameters:
491: .  rg - the region context

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

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

502:    Level: intermediate

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


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

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

533:    Logically Collective on RG

535:    Input Parameters:
536: +  rg  - the region context
537: -  flg - the boolean flag

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

542:    Level: intermediate

544: .seealso: RGGetComplement()
545: @*/
546: PetscErrorCode RGSetComplement(RG rg,PetscBool flg)
547: {
551:   rg->complement = flg;
552:   return(0);
553: }

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

559:    Not Collective

561:    Input Parameter:
562: .  rg - the region context

564:    Output Parameter:
565: .  flg - the flag

567:    Level: intermediate

569: .seealso: RGSetComplement()
570: @*/
571: PetscErrorCode RGGetComplement(RG rg,PetscBool *flg)
572: {
576:   *flg = rg->complement;
577:   return(0);
578: }

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

584:    Logically Collective on RG

586:    Input Parameters:
587: +  rg      - the region context
588: -  sfactor - the scaling factor

590:    Options Database Key:
591: .  -rg_scale <real> - Sets the scaling factor

593:    Level: advanced

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

610: /*@
611:    RGGetScale - Gets the scaling factor.

613:    Not Collective

615:    Input Parameter:
616: .  rg - the region context

618:    Output Parameter:
619: .  flg - the flag

621:    Level: advanced

623: .seealso: RGSetScale()
624: @*/
625: PetscErrorCode RGGetScale(RG rg,PetscReal *sfactor)
626: {
630:   *sfactor = rg->sfactor;
631:   return(0);
632: }

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

638:    Logically Collective on RG

640:    Input Parameters:
641: +  rg      - the region context
642: -  sfactor - the scaling factor

644:    Notes:
645:    The current implementation does not allow pushing several scaling factors.

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

650:    Level: developer

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

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

669:    Not Collective

671:    Input Parameter:
672: .  rg - the region context

674:    Level: developer

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

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

691:    Collective on RG

693:    Input Parameter:
694: .  rg - the region context

696:    Level: beginner

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

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

713: /*@C
714:    RGRegister - Adds a region to the RG package.

716:    Not collective

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

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

725:    Level: advanced

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

734:   RGInitializePackage();
735:   PetscFunctionListAdd(&RGList,name,function);
736:   return(0);
737: }