Actual source code: slepcrg.h

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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:    User interface for the region object in SLEPc
 12: */

 14: #pragma once

 16: #include <slepcsys.h>
 17: #include <slepcrgtypes.h>

 19: /* SUBMANSEC = RG */

 21: SLEPC_EXTERN PetscErrorCode RGInitializePackage(void);
 22: SLEPC_EXTERN PetscErrorCode RGFinalizePackage(void);

 24: /*J
 25:    RGType - String with the name of the region type.

 27:    Level: beginner

 29: .seealso: [](sec:rg), `RG`, `RGSetType()`
 30: J*/
 31: typedef const char *RGType;
 32: #define RGINTERVAL  "interval"
 33: #define RGPOLYGON   "polygon"
 34: #define RGELLIPSE   "ellipse"
 35: #define RGRING      "ring"

 37: /* Logging support */
 38: SLEPC_EXTERN PetscClassId RG_CLASSID;

 40: SLEPC_EXTERN PetscErrorCode RGCreate(MPI_Comm,RG*);
 41: SLEPC_EXTERN PetscErrorCode RGSetType(RG,RGType);
 42: SLEPC_EXTERN PetscErrorCode RGGetType(RG,RGType*);
 43: SLEPC_EXTERN PetscErrorCode RGSetOptionsPrefix(RG,const char[]);
 44: SLEPC_EXTERN PetscErrorCode RGAppendOptionsPrefix(RG,const char[]);
 45: SLEPC_EXTERN PetscErrorCode RGGetOptionsPrefix(RG,const char*[]);
 46: SLEPC_EXTERN PetscErrorCode RGSetFromOptions(RG);
 47: SLEPC_EXTERN PetscErrorCode RGView(RG,PetscViewer);
 48: SLEPC_EXTERN PetscErrorCode RGViewFromOptions(RG,PetscObject,const char[]);
 49: SLEPC_EXTERN PetscErrorCode RGDestroy(RG*);

 51: /*E
 52:    RGQuadRule - Determines how to compute the quadrature parameters.

 54:    Values:
 55: +  `RG_QUADRULE_TRAPEZOIDAL` - use the trapezoidal rule
 56: -  `RG_QUADRULE_CHEBYSHEV`   - use Chebyshev nodes and weights

 58:    Level: advanced

 60: .seealso: [](sec:rg), `RGComputeQuadrature()`
 61: E*/
 62: typedef enum { RG_QUADRULE_TRAPEZOIDAL = 1,
 63:                RG_QUADRULE_CHEBYSHEV   = 2 } RGQuadRule;

 65: /*MC
 66:    RG_QUADRULE_TRAPEZOIDAL - In contour integral methods, use the trapezoidal rule.

 68:    Level: advanced

 70: .seealso: [](sec:rg), `RGQuadRule`, `RGComputeQuadrature()`, `RG_QUADRULE_CHEBYSHEV`
 71: M*/

 73: /*MC
 74:    RG_QUADRULE_CHEBYSHEV - In contour integral methods, use Chebyshev nodes and weights.

 76:    Level: advanced

 78: .seealso: [](sec:rg), `RGQuadRule`, `RGComputeQuadrature()`, `RG_QUADRULE_TRAPEZOIDAL`
 79: M*/

 81: SLEPC_EXTERN PetscErrorCode RGIsTrivial(RG,PetscBool*);
 82: SLEPC_EXTERN PetscErrorCode RGSetComplement(RG,PetscBool);
 83: SLEPC_EXTERN PetscErrorCode RGGetComplement(RG,PetscBool*);
 84: SLEPC_EXTERN PetscErrorCode RGSetScale(RG,PetscReal);
 85: SLEPC_EXTERN PetscErrorCode RGGetScale(RG,PetscReal*);
 86: SLEPC_EXTERN PetscErrorCode RGPushScale(RG,PetscReal);
 87: SLEPC_EXTERN PetscErrorCode RGPopScale(RG);
 88: SLEPC_EXTERN PetscErrorCode RGCheckInside(RG,PetscInt,PetscScalar[],PetscScalar[],PetscInt[]);
 89: SLEPC_EXTERN PetscErrorCode RGIsAxisymmetric(RG,PetscBool,PetscBool*);
 90: SLEPC_EXTERN PetscErrorCode RGCanUseConjugates(RG,PetscBool,PetscBool*);
 91: SLEPC_EXTERN PetscErrorCode RGComputeContour(RG,PetscInt,PetscScalar[],PetscScalar[]);
 92: SLEPC_EXTERN PetscErrorCode RGComputeBoundingBox(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
 93: SLEPC_EXTERN PetscErrorCode RGComputeQuadrature(RG,RGQuadRule,PetscInt,PetscScalar[],PetscScalar[],PetscScalar[]);

 95: SLEPC_EXTERN PetscFunctionList RGList;
 96: SLEPC_EXTERN PetscErrorCode RGRegister(const char[],PetscErrorCode(*)(RG));

 98: /* --------- options specific to particular regions -------- */

100: SLEPC_EXTERN PetscErrorCode RGEllipseSetParameters(RG,PetscScalar,PetscReal,PetscReal);
101: SLEPC_EXTERN PetscErrorCode RGEllipseGetParameters(RG,PetscScalar*,PetscReal*,PetscReal*);

103: SLEPC_EXTERN PetscErrorCode RGIntervalSetEndpoints(RG,PetscReal,PetscReal,PetscReal,PetscReal);
104: SLEPC_EXTERN PetscErrorCode RGIntervalGetEndpoints(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);

106: SLEPC_EXTERN PetscErrorCode RGPolygonSetVertices(RG,PetscInt,PetscScalar[],PetscScalar[]);
107: SLEPC_EXTERN PetscErrorCode RGPolygonGetVertices(RG,PetscInt*,PetscScalar*[],PetscScalar*[]);

109: SLEPC_EXTERN PetscErrorCode RGRingSetParameters(RG,PetscScalar,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
110: SLEPC_EXTERN PetscErrorCode RGRingGetParameters(RG,PetscScalar*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*);