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: Basic RG routines
12: */
14: #include <slepc/private/rgimpl.h> 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;
53: PetscClassId classids[1];
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 */
64: classids[0] = RG_CLASSID;
65: PetscInfoProcessClass("rg",1,&classids[0]);
66: /* Process summary exclusions */
67: PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);
68: if (opt) {
69: PetscStrInList("rg",logList,',',&pkg);
70: if (pkg) { PetscLogEventDeactivateClass(RG_CLASSID); }
71: }
72: /* Register package finalizer */
73: PetscRegisterFinalize(RGFinalizePackage);
74: return(0);
75: }
77: /*@
78: RGCreate - Creates an RG context.
80: Collective
82: Input Parameter:
83: . comm - MPI communicator
85: Output Parameter:
86: . newrg - location to put the RG context
88: Level: beginner
90: .seealso: RGDestroy(), RG 91: @*/
92: PetscErrorCode RGCreate(MPI_Comm comm,RG *newrg) 93: {
94: RG rg;
99: *newrg = 0;
100: RGInitializePackage();
101: SlepcHeaderCreate(rg,RG_CLASSID,"RG","Region","RG",comm,RGDestroy,RGView);
102: rg->complement = PETSC_FALSE;
103: rg->sfactor = 1.0;
104: rg->osfactor = 0.0;
105: rg->data = NULL;
107: *newrg = rg;
108: return(0);
109: }
111: /*@C
112: RGSetOptionsPrefix - Sets the prefix used for searching for all
113: RG options in the database.
115: Logically Collective on rg
117: Input Parameters:
118: + rg - the region context
119: - prefix - the prefix string to prepend to all RG option requests
121: Notes:
122: A hyphen (-) must NOT be given at the beginning of the prefix name.
123: The first character of all runtime options is AUTOMATICALLY the
124: hyphen.
126: Level: advanced
128: .seealso: RGAppendOptionsPrefix()
129: @*/
130: PetscErrorCode RGSetOptionsPrefix(RG rg,const char *prefix)131: {
136: PetscObjectSetOptionsPrefix((PetscObject)rg,prefix);
137: return(0);
138: }
140: /*@C
141: RGAppendOptionsPrefix - Appends to the prefix used for searching for all
142: RG options in the database.
144: Logically Collective on rg
146: Input Parameters:
147: + rg - the region context
148: - prefix - the prefix string to prepend to all RG option requests
150: Notes:
151: A hyphen (-) must NOT be given at the beginning of the prefix name.
152: The first character of all runtime options is AUTOMATICALLY the hyphen.
154: Level: advanced
156: .seealso: RGSetOptionsPrefix()
157: @*/
158: PetscErrorCode RGAppendOptionsPrefix(RG rg,const char *prefix)159: {
164: PetscObjectAppendOptionsPrefix((PetscObject)rg,prefix);
165: return(0);
166: }
168: /*@C
169: RGGetOptionsPrefix - Gets the prefix used for searching for all
170: RG options in the database.
172: Not Collective
174: Input Parameters:
175: . rg - the region context
177: Output Parameters:
178: . prefix - pointer to the prefix string used is returned
180: Note:
181: On the Fortran side, the user should pass in a string 'prefix' of
182: sufficient length to hold the prefix.
184: Level: advanced
186: .seealso: RGSetOptionsPrefix(), RGAppendOptionsPrefix()
187: @*/
188: PetscErrorCode RGGetOptionsPrefix(RG rg,const char *prefix[])189: {
195: PetscObjectGetOptionsPrefix((PetscObject)rg,prefix);
196: return(0);
197: }
199: /*@C
200: RGSetType - Selects the type for the RG object.
202: Logically Collective on rg
204: Input Parameters:
205: + rg - the region context
206: - type - a known type
208: Level: intermediate
210: .seealso: RGGetType()
211: @*/
212: PetscErrorCode RGSetType(RG rg,RGType type)213: {
214: PetscErrorCode ierr,(*r)(RG);
215: PetscBool match;
221: PetscObjectTypeCompare((PetscObject)rg,type,&match);
222: if (match) return(0);
224: PetscFunctionListFind(RGList,type,&r);
225: if (!r) SETERRQ1(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested RG type %s",type);
227: if (rg->ops->destroy) { (*rg->ops->destroy)(rg); }
228: PetscMemzero(rg->ops,sizeof(struct _RGOps));
230: PetscObjectChangeTypeName((PetscObject)rg,type);
231: (*r)(rg);
232: return(0);
233: }
235: /*@C
236: RGGetType - Gets the RG type name (as a string) from the RG context.
238: Not Collective
240: Input Parameter:
241: . rg - the region context
243: Output Parameter:
244: . name - name of the region
246: Level: intermediate
248: .seealso: RGSetType()
249: @*/
250: PetscErrorCode RGGetType(RG rg,RGType *type)251: {
255: *type = ((PetscObject)rg)->type_name;
256: return(0);
257: }
259: /*@
260: RGSetFromOptions - Sets RG options from the options database.
262: Collective on rg
264: Input Parameters:
265: . rg - the region context
267: Notes:
268: To see all options, run your program with the -help option.
270: Level: beginner
271: @*/
272: PetscErrorCode RGSetFromOptions(RG rg)273: {
275: char type[256];
276: PetscBool flg;
277: PetscReal sfactor;
281: RGRegisterAll();
282: PetscObjectOptionsBegin((PetscObject)rg);
283: PetscOptionsFList("-rg_type","Region type","RGSetType",RGList,(char*)(((PetscObject)rg)->type_name?((PetscObject)rg)->type_name:RGINTERVAL),type,sizeof(type),&flg);
284: if (flg) {
285: RGSetType(rg,type);
286: } else if (!((PetscObject)rg)->type_name) {
287: RGSetType(rg,RGINTERVAL);
288: }
290: PetscOptionsBool("-rg_complement","Whether region is complemented or not","RGSetComplement",rg->complement,&rg->complement,NULL);
292: PetscOptionsReal("-rg_scale","Scaling factor","RGSetScale",1.0,&sfactor,&flg);
293: if (flg) { RGSetScale(rg,sfactor); }
295: if (rg->ops->setfromoptions) {
296: (*rg->ops->setfromoptions)(PetscOptionsObject,rg);
297: }
298: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)rg);
299: PetscOptionsEnd();
300: return(0);
301: }
303: /*@C
304: RGView - Prints the RG data structure.
306: Collective on rg
308: Input Parameters:
309: + rg - the region context
310: - viewer - optional visualization context
312: Note:
313: The available visualization contexts include
314: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
315: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
316: output where only the first processor opens
317: the file. All other processors send their
318: data to the first processor to print.
320: The user can open an alternative visualization context with
321: PetscViewerASCIIOpen() - output to a specified file.
323: Level: beginner
324: @*/
325: PetscErrorCode RGView(RG rg,PetscViewer viewer)326: {
327: PetscBool isdraw,isascii;
332: if (!viewer) {
333: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)rg),&viewer);
334: }
337: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
338: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
339: if (isascii) {
340: PetscObjectPrintClassNamePrefixType((PetscObject)rg,viewer);
341: if (rg->ops->view) {
342: PetscViewerASCIIPushTab(viewer);
343: (*rg->ops->view)(rg,viewer);
344: PetscViewerASCIIPopTab(viewer);
345: }
346: if (rg->complement) {
347: PetscViewerASCIIPrintf(viewer," selected region is the complement of the specified one\n");
348: }
349: if (rg->sfactor!=1.0) {
350: PetscViewerASCIIPrintf(viewer," scaling factor = %g\n",(double)rg->sfactor);
351: }
352: } else if (isdraw) {
353: if (rg->ops->view) { (*rg->ops->view)(rg,viewer); }
354: }
355: return(0);
356: }
358: /*@C
359: RGViewFromOptions - View from options
361: Collective on RG363: Input Parameters:
364: + rg - the region context
365: . obj - optional object
366: - name - command line option
368: Level: intermediate
370: .seealso: RGView(), RGCreate()
371: @*/
372: PetscErrorCode RGViewFromOptions(RG rg,PetscObject obj,const char name[])373: {
378: PetscObjectViewFromOptions((PetscObject)rg,obj,name);
379: return(0);
380: }
382: /*@
383: RGIsTrivial - Whether it is the trivial region (whole complex plane).
385: Not Collective
387: Input Parameter:
388: . rg - the region context
390: Output Parameter:
391: . trivial - true if the region is equal to the whole complex plane, e.g.,
392: an interval region with all four endpoints unbounded or an
393: ellipse with infinite radius.
395: Level: beginner
396: @*/
397: PetscErrorCode RGIsTrivial(RG rg,PetscBool *trivial)398: {
405: if (rg->ops->istrivial) {
406: (*rg->ops->istrivial)(rg,trivial);
407: } else *trivial = PETSC_FALSE;
408: return(0);
409: }
411: /*@
412: RGCheckInside - Determines if a set of given points are inside the region or not.
414: Not Collective
416: Input Parameters:
417: + rg - the region context
418: . n - number of points to check
419: . ar - array of real parts
420: - ai - array of imaginary parts
422: Output Parameter:
423: . inside - array of results (1=inside, 0=on the contour, -1=outside)
425: Note:
426: The point a is expressed as a couple of PetscScalar variables ar,ai.
427: If built with complex scalars, the point is supposed to be stored in ar,
428: otherwise ar,ai contain the real and imaginary parts, respectively.
430: If a scaling factor was set, the points are scaled before checking.
432: Level: intermediate
434: .seealso: RGSetScale(), RGSetComplement()
435: @*/
436: PetscErrorCode RGCheckInside(RG rg,PetscInt n,PetscScalar *ar,PetscScalar *ai,PetscInt *inside)437: {
439: PetscReal px,py;
440: PetscInt i;
446: #if !defined(PETSC_USE_COMPLEX)
448: #endif
451: for (i=0;i<n;i++) {
452: #if defined(PETSC_USE_COMPLEX)
453: px = PetscRealPart(ar[i]);
454: py = PetscImaginaryPart(ar[i]);
455: #else
456: px = ar[i];
457: py = ai[i];
458: #endif
459: if (rg->sfactor != 1.0) {
460: px /= rg->sfactor;
461: py /= rg->sfactor;
462: }
463: (*rg->ops->checkinside)(rg,px,py,inside+i);
464: if (rg->complement) inside[i] = -inside[i];
465: }
466: return(0);
467: }
469: /*@
470: RGIsAxisymmetric - Determines if the region is symmetric with respect
471: to the real or imaginary axis.
473: Not Collective
475: Input Parameters:
476: + rg - the region context
477: - vertical - true if symmetry must be checked against the vertical axis
479: Output Parameter:
480: . symm - true if the region is axisymmetric
482: Note:
483: If the vertical argument is true, symmetry is checked with respect to
484: the vertical axis, otherwise with respect to the horizontal axis.
486: Level: intermediate
487: @*/
488: PetscErrorCode RGIsAxisymmetric(RG rg,PetscBool vertical,PetscBool *symm)489: {
497: if (rg->ops->isaxisymmetric) {
498: (*rg->ops->isaxisymmetric)(rg,vertical,symm);
499: } else *symm = PETSC_FALSE;
500: return(0);
501: }
503: /*@
504: RGCanUseConjugates - Used in contour integral methods to determine whether
505: half of integration points can be avoided (use their conjugates).
507: Not Collective
509: Input Parameters:
510: + rg - the region context
511: - realmats - true if the problem matrices are real
513: Output Parameter:
514: . useconj - whether it is possible to use conjugates
516: Notes:
517: If some integration points are the conjugates of other points, then the
518: associated computational cost can be saved. This depends on the problem
519: matrices being real and also the region being symmetric with respect to
520: the horizontal axis. The result is false if using real arithmetic or
521: in the case of a flat region (height equal to zero).
523: Level: developer
524: @*/
525: PetscErrorCode RGCanUseConjugates(RG rg,PetscBool realmats,PetscBool *useconj)526: {
527: #if defined(PETSC_USE_COMPLEX)
529: PetscReal c,d;
530: PetscBool isaxisymm;
531: #endif
537: *useconj = PETSC_FALSE;
538: #if defined(PETSC_USE_COMPLEX)
539: if (realmats) {
540: RGIsAxisymmetric(rg,PETSC_FALSE,&isaxisymm);
541: if (isaxisymm) {
542: RGComputeBoundingBox(rg,NULL,NULL,&c,&d);
543: if (c!=d) *useconj = PETSC_TRUE;
544: }
545: }
546: #endif
547: return(0);
548: }
550: /*@
551: RGComputeContour - Computes the coordinates of several points lying on the
552: contour of the region.
554: Not Collective
556: Input Parameters:
557: + rg - the region context
558: - n - number of points to compute
560: Output Parameters:
561: + cr - location to store real parts
562: - ci - location to store imaginary parts
564: Notes:
565: In real scalars, either cr or ci can be NULL (but not both). In complex
566: scalars, the coordinates are stored in cr, which cannot be NULL (ci is
567: not referenced).
569: Level: intermediate
570: @*/
571: PetscErrorCode RGComputeContour(RG rg,PetscInt n,PetscScalar cr[],PetscScalar ci[])572: {
573: PetscInt i;
579: #if defined(PETSC_USE_COMPLEX)
581: #else
582: if (!cr && !ci) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"cr and ci cannot be NULL simultaneously");
583: #endif
584: if (rg->complement) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Cannot compute contour of region with complement flag set");
585: (*rg->ops->computecontour)(rg,n,cr,ci);
586: for (i=0;i<n;i++) {
587: if (cr) cr[i] *= rg->sfactor;
588: if (ci) ci[i] *= rg->sfactor;
589: }
590: return(0);
591: }
593: /*@
594: RGComputeBoundingBox - Determines the endpoints of a rectangle in the complex plane that
595: contains the region.
597: Not Collective
599: Input Parameter:
600: . rg - the region context
602: Output Parameters:
603: + a - left endpoint of the bounding box in the real axis
604: . b - right endpoint of the bounding box in the real axis
605: . c - bottom endpoint of the bounding box in the imaginary axis
606: - d - top endpoint of the bounding box in the imaginary axis
608: Notes:
609: The bounding box is defined as [a,b]x[c,d]. In regions that are not bounded (e.g. an
610: open interval) or with the complement flag set, it makes no sense to compute a bounding
611: box, so the return values are infinite.
613: Level: intermediate
615: .seealso: RGSetComplement()
616: @*/
617: PetscErrorCode RGComputeBoundingBox(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)618: {
625: if (rg->complement) { /* cannot compute bounding box */
626: if (a) *a = -PETSC_MAX_REAL;
627: if (b) *b = PETSC_MAX_REAL;
628: if (c) *c = -PETSC_MAX_REAL;
629: if (d) *d = PETSC_MAX_REAL;
630: } else {
631: (*rg->ops->computebbox)(rg,a,b,c,d);
632: if (a && *a!=-PETSC_MAX_REAL) *a *= rg->sfactor;
633: if (b && *b!= PETSC_MAX_REAL) *b *= rg->sfactor;
634: if (c && *c!=-PETSC_MAX_REAL) *c *= rg->sfactor;
635: if (d && *d!= PETSC_MAX_REAL) *d *= rg->sfactor;
636: }
637: return(0);
638: }
640: /*@
641: RGComputeQuadrature - Computes the values of the parameters used in a
642: quadrature rule for a contour integral around the boundary of the region.
644: Not Collective
646: Input Parameters:
647: + rg - the region context
648: . quad - the type of quadrature
649: - n - number of quadrature points to compute
651: Output Parameters:
652: + z - quadrature points
653: . zn - normalized quadrature points
654: - w - quadrature weights
656: Notes:
657: In complex scalars, the values returned in z are often the same as those
658: computed by RGComputeContour(), but this is not the case in real scalars
659: where all output arguments are real.
661: The computed values change for different quadrature rules (trapezoidal
662: or Chebyshev).
664: Level: intermediate
666: .seealso: RGComputeContour()
667: @*/
668: PetscErrorCode RGComputeQuadrature(RG rg,RGQuadRule quad,PetscInt n,PetscScalar z[],PetscScalar zn[],PetscScalar w[])669: {
679: RGComputeContour(rg,n,z,NULL);
680: if (!rg->ops->computequadrature) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Quadrature rules not implemented for this region type");
681: (*rg->ops->computequadrature)(rg,quad,n,z,zn,w);
682: return(0);
683: }
685: /*@
686: RGSetComplement - Sets a flag to indicate that the region is the complement
687: of the specified one.
689: Logically Collective on rg
691: Input Parameters:
692: + rg - the region context
693: - flg - the boolean flag
695: Options Database Key:
696: . -rg_complement <bool> - Activate/deactivate the complementation of the region
698: Level: intermediate
700: .seealso: RGGetComplement()
701: @*/
702: PetscErrorCode RGSetComplement(RG rg,PetscBool flg)703: {
707: rg->complement = flg;
708: return(0);
709: }
711: /*@
712: RGGetComplement - Gets a flag that that indicates whether the region
713: is complemented or not.
715: Not Collective
717: Input Parameter:
718: . rg - the region context
720: Output Parameter:
721: . flg - the flag
723: Level: intermediate
725: .seealso: RGSetComplement()
726: @*/
727: PetscErrorCode RGGetComplement(RG rg,PetscBool *flg)728: {
732: *flg = rg->complement;
733: return(0);
734: }
736: /*@
737: RGSetScale - Sets the scaling factor to be used when checking that a
738: point is inside the region and when computing the contour.
740: Logically Collective on rg
742: Input Parameters:
743: + rg - the region context
744: - sfactor - the scaling factor
746: Options Database Key:
747: . -rg_scale <real> - Sets the scaling factor
749: Level: advanced
751: .seealso: RGGetScale(), RGCheckInside()
752: @*/
753: PetscErrorCode RGSetScale(RG rg,PetscReal sfactor)754: {
758: if (sfactor == PETSC_DEFAULT || sfactor == PETSC_DECIDE) rg->sfactor = 1.0;
759: else {
760: if (sfactor<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of scaling factor. Must be > 0");
761: rg->sfactor = sfactor;
762: }
763: return(0);
764: }
766: /*@
767: RGGetScale - Gets the scaling factor.
769: Not Collective
771: Input Parameter:
772: . rg - the region context
774: Output Parameter:
775: . flg - the flag
777: Level: advanced
779: .seealso: RGSetScale()
780: @*/
781: PetscErrorCode RGGetScale(RG rg,PetscReal *sfactor)782: {
786: *sfactor = rg->sfactor;
787: return(0);
788: }
790: /*@
791: RGPushScale - Sets an additional scaling factor, that will multiply the
792: user-defined scaling factor.
794: Logically Collective on rg
796: Input Parameters:
797: + rg - the region context
798: - sfactor - the scaling factor
800: Notes:
801: The current implementation does not allow pushing several scaling factors.
803: This is intended for internal use, for instance in polynomial eigensolvers
804: that use parameter scaling.
806: Level: developer
808: .seealso: RGPopScale(), RGSetScale()
809: @*/
810: PetscErrorCode RGPushScale(RG rg,PetscReal sfactor)811: {
815: if (sfactor<=0.0) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of scaling factor. Must be > 0");
816: if (rg->osfactor) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Current implementation does not allow pushing several scaling factors");
817: rg->osfactor = rg->sfactor;
818: rg->sfactor *= sfactor;
819: return(0);
820: }
822: /*@
823: RGPopScale - Pops the scaling factor set with RGPushScale().
825: Not Collective
827: Input Parameter:
828: . rg - the region context
830: Level: developer
832: .seealso: RGPushScale()
833: @*/
834: PetscErrorCode RGPopScale(RG rg)835: {
838: if (!rg->osfactor) SETERRQ(PetscObjectComm((PetscObject)rg),PETSC_ERR_ORDER,"Must call RGPushScale first");
839: rg->sfactor = rg->osfactor;
840: rg->osfactor = 0.0;
841: return(0);
842: }
844: /*@C
845: RGDestroy - Destroys RG context that was created with RGCreate().
847: Collective on rg
849: Input Parameter:
850: . rg - the region context
852: Level: beginner
854: .seealso: RGCreate()
855: @*/
856: PetscErrorCode RGDestroy(RG *rg)857: {
861: if (!*rg) return(0);
863: if (--((PetscObject)(*rg))->refct > 0) { *rg = 0; return(0); }
864: if ((*rg)->ops->destroy) { (*(*rg)->ops->destroy)(*rg); }
865: PetscHeaderDestroy(rg);
866: return(0);
867: }
869: /*@C
870: RGRegister - Adds a region to the RG package.
872: Not collective
874: Input Parameters:
875: + name - name of a new user-defined RG876: - function - routine to create context
878: Notes:
879: RGRegister() may be called multiple times to add several user-defined regions.
881: Level: advanced
883: .seealso: RGRegisterAll()
884: @*/
885: PetscErrorCode RGRegister(const char *name,PetscErrorCode (*function)(RG))886: {
890: RGInitializePackage();
891: PetscFunctionListAdd(&RGList,name,function);
892: return(0);
893: }