Actual source code: slepcrg.h
slepc-3.17.1 2022-04-11
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: #if !defined(SLEPCRG_H)
15: #define SLEPCRG_H
16: #include <slepcsys.h>
17: #include <slepcrgtypes.h>
19: SLEPC_EXTERN PetscErrorCode RGInitializePackage(void);
21: /*J
22: RGType - String with the name of the region.
24: Level: beginner
26: .seealso: RGSetType(), RG
27: J*/
28: typedef const char* RGType;
29: #define RGINTERVAL "interval"
30: #define RGPOLYGON "polygon"
31: #define RGELLIPSE "ellipse"
32: #define RGRING "ring"
34: /* Logging support */
35: SLEPC_EXTERN PetscClassId RG_CLASSID;
37: SLEPC_EXTERN PetscErrorCode RGCreate(MPI_Comm,RG*);
38: SLEPC_EXTERN PetscErrorCode RGSetType(RG,RGType);
39: SLEPC_EXTERN PetscErrorCode RGGetType(RG,RGType*);
40: SLEPC_EXTERN PetscErrorCode RGSetOptionsPrefix(RG,const char *);
41: SLEPC_EXTERN PetscErrorCode RGAppendOptionsPrefix(RG,const char *);
42: SLEPC_EXTERN PetscErrorCode RGGetOptionsPrefix(RG,const char *[]);
43: SLEPC_EXTERN PetscErrorCode RGSetFromOptions(RG);
44: SLEPC_EXTERN PetscErrorCode RGView(RG,PetscViewer);
45: SLEPC_EXTERN PetscErrorCode RGViewFromOptions(RG,PetscObject,const char[]);
46: SLEPC_EXTERN PetscErrorCode RGDestroy(RG*);
48: /*E
49: RGQuadRule - determines how to compute the quadrature parameters
51: Level: advanced
53: .seealso: RGComputeQuadrature()
54: E*/
55: typedef enum { RG_QUADRULE_TRAPEZOIDAL=1,
56: RG_QUADRULE_CHEBYSHEV } RGQuadRule;
58: SLEPC_EXTERN PetscErrorCode RGIsTrivial(RG,PetscBool*);
59: SLEPC_EXTERN PetscErrorCode RGSetComplement(RG,PetscBool);
60: SLEPC_EXTERN PetscErrorCode RGGetComplement(RG,PetscBool*);
61: SLEPC_EXTERN PetscErrorCode RGSetScale(RG,PetscReal);
62: SLEPC_EXTERN PetscErrorCode RGGetScale(RG,PetscReal*);
63: SLEPC_EXTERN PetscErrorCode RGPushScale(RG,PetscReal);
64: SLEPC_EXTERN PetscErrorCode RGPopScale(RG);
65: SLEPC_EXTERN PetscErrorCode RGCheckInside(RG,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);
66: SLEPC_EXTERN PetscErrorCode RGIsAxisymmetric(RG,PetscBool,PetscBool*);
67: SLEPC_EXTERN PetscErrorCode RGCanUseConjugates(RG,PetscBool,PetscBool*);
68: SLEPC_EXTERN PetscErrorCode RGComputeContour(RG,PetscInt,PetscScalar*,PetscScalar*);
69: SLEPC_EXTERN PetscErrorCode RGComputeBoundingBox(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
70: SLEPC_EXTERN PetscErrorCode RGComputeQuadrature(RG,RGQuadRule,PetscInt,PetscScalar*,PetscScalar*,PetscScalar*);
72: SLEPC_EXTERN PetscFunctionList RGList;
73: SLEPC_EXTERN PetscErrorCode RGRegister(const char[],PetscErrorCode(*)(RG));
75: /* --------- options specific to particular regions -------- */
77: SLEPC_EXTERN PetscErrorCode RGEllipseSetParameters(RG,PetscScalar,PetscReal,PetscReal);
78: SLEPC_EXTERN PetscErrorCode RGEllipseGetParameters(RG,PetscScalar*,PetscReal*,PetscReal*);
80: SLEPC_EXTERN PetscErrorCode RGIntervalSetEndpoints(RG,PetscReal,PetscReal,PetscReal,PetscReal);
81: SLEPC_EXTERN PetscErrorCode RGIntervalGetEndpoints(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
83: SLEPC_EXTERN PetscErrorCode RGPolygonSetVertices(RG,PetscInt,PetscScalar*,PetscScalar*);
84: SLEPC_EXTERN PetscErrorCode RGPolygonGetVertices(RG,PetscInt*,PetscScalar**,PetscScalar**);
86: SLEPC_EXTERN PetscErrorCode RGRingSetParameters(RG,PetscScalar,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
87: SLEPC_EXTERN PetscErrorCode RGRingGetParameters(RG,PetscScalar*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
89: #endif