Actual source code: rgbasic.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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:   char             *className;
 53:   PetscBool        opt;
 54:   PetscErrorCode   ierr;

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

 83: /*@
 84:    RGCreate - Creates an RG context.

 86:    Collective on MPI_Comm

 88:    Input Parameter:
 89: .  comm - MPI communicator

 91:    Output Parameter:
 92: .  newrg - location to put the RG context

 94:    Level: beginner

 96: .seealso: RGDestroy(), RG
 97: @*/
 98: PetscErrorCode RGCreate(MPI_Comm comm,RG *newrg)
 99: {
100:   RG             rg;

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

113:   *newrg = rg;
114:   return(0);
115: }

117: /*@C
118:    RGSetOptionsPrefix - Sets the prefix used for searching for all
119:    RG options in the database.

121:    Logically Collective on RG

123:    Input Parameters:
124: +  rg     - the region context
125: -  prefix - the prefix string to prepend to all RG option requests

127:    Notes:
128:    A hyphen (-) must NOT be given at the beginning of the prefix name.
129:    The first character of all runtime options is AUTOMATICALLY the
130:    hyphen.

132:    Level: advanced

134: .seealso: RGAppendOptionsPrefix()
135: @*/
136: PetscErrorCode RGSetOptionsPrefix(RG rg,const char *prefix)
137: {

142:   PetscObjectSetOptionsPrefix((PetscObject)rg,prefix);
143:   return(0);
144: }

146: /*@C
147:    RGAppendOptionsPrefix - Appends to the prefix used for searching for all
148:    RG options in the database.

150:    Logically Collective on RG

152:    Input Parameters:
153: +  rg     - the region context
154: -  prefix - the prefix string to prepend to all RG option requests

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

160:    Level: advanced

162: .seealso: RGSetOptionsPrefix()
163: @*/
164: PetscErrorCode RGAppendOptionsPrefix(RG rg,const char *prefix)
165: {

170:   PetscObjectAppendOptionsPrefix((PetscObject)rg,prefix);
171:   return(0);
172: }

174: /*@C
175:    RGGetOptionsPrefix - Gets the prefix used for searching for all
176:    RG options in the database.

178:    Not Collective

180:    Input Parameters:
181: .  rg - the region context

183:    Output Parameters:
184: .  prefix - pointer to the prefix string used is returned

186:    Note:
187:    On the Fortran side, the user should pass in a string 'prefix' of
188:    sufficient length to hold the prefix.

190:    Level: advanced

192: .seealso: RGSetOptionsPrefix(), RGAppendOptionsPrefix()
193: @*/
194: PetscErrorCode RGGetOptionsPrefix(RG rg,const char *prefix[])
195: {

201:   PetscObjectGetOptionsPrefix((PetscObject)rg,prefix);
202:   return(0);
203: }

205: /*@C
206:    RGSetType - Selects the type for the RG object.

208:    Logically Collective on RG

210:    Input Parameter:
211: +  rg   - the region context
212: -  type - a known type

214:    Level: intermediate

216: .seealso: RGGetType()
217: @*/
218: PetscErrorCode RGSetType(RG rg,RGType type)
219: {
220:   PetscErrorCode ierr,(*r)(RG);
221:   PetscBool      match;


227:   PetscObjectTypeCompare((PetscObject)rg,type,&match);
228:   if (match) return(0);

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

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

236:   PetscObjectChangeTypeName((PetscObject)rg,type);
237:   (*r)(rg);
238:   return(0);
239: }

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

244:    Not Collective

246:    Input Parameter:
247: .  rg - the region context

249:    Output Parameter:
250: .  name - name of the region

252:    Level: intermediate

254: .seealso: RGSetType()
255: @*/
256: PetscErrorCode RGGetType(RG rg,RGType *type)
257: {
261:   *type = ((PetscObject)rg)->type_name;
262:   return(0);
263: }

265: /*@
266:    RGSetFromOptions - Sets RG options from the options database.

268:    Collective on RG

270:    Input Parameters:
271: .  rg - the region context

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

276:    Level: beginner
277: @*/
278: PetscErrorCode RGSetFromOptions(RG rg)
279: {
281:   char           type[256];
282:   PetscBool      flg;
283:   PetscReal      sfactor;

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

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

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

301:     if (rg->ops->setfromoptions) {
302:       (*rg->ops->setfromoptions)(PetscOptionsObject,rg);
303:     }
304:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)rg);
305:   PetscOptionsEnd();
306:   return(0);
307: }

309: /*@C
310:    RGView - Prints the RG data structure.

312:    Collective on RG

314:    Input Parameters:
315: +  rg - the region context
316: -  viewer - optional visualization context

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

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

329:    Level: beginner
330: @*/
331: PetscErrorCode RGView(RG rg,PetscViewer viewer)
332: {
333:   PetscBool      isascii;

338:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)rg));
341:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
342:   if (isascii) {
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:   }
356:   return(0);
357: }

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

362:    Not Collective

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

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

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

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

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

391:    Not Collective

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

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

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

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

409:    Level: intermediate

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

423: #if !defined(PETSC_USE_COMPLEX)
425: #endif

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

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

450:    Not Collective

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

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

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

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

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

487:    Not Collective

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

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

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

501:    Level: intermediate

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


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

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

532:    Logically Collective on RG

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

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

541:    Level: intermediate

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

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

558:    Not Collective

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

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

566:    Level: intermediate

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

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

583:    Logically Collective on RG

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

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

592:    Level: advanced

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

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

612:    Not Collective

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

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

620:    Level: advanced

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

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

637:    Logically Collective on RG

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

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

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

649:    Level: developer

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

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

668:    Not Collective

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

673:    Level: developer

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

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

690:    Collective on RG

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

695:    Level: beginner

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

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

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

715:    Not collective

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

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

724:    Level: advanced

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

733:   PetscFunctionListAdd(&RGList,name,function);
734:   return(0);
735: }