Actual source code: slepcrg.h
slepc-3.22.1 2024-10-28
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.
27: Level: beginner
29: .seealso: RGSetType(), RG
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: Level: advanced
56: .seealso: RGComputeQuadrature()
57: E*/
58: typedef enum { RG_QUADRULE_TRAPEZOIDAL=1,
59: RG_QUADRULE_CHEBYSHEV } RGQuadRule;
61: SLEPC_EXTERN PetscErrorCode RGIsTrivial(RG,PetscBool*);
62: SLEPC_EXTERN PetscErrorCode RGSetComplement(RG,PetscBool);
63: SLEPC_EXTERN PetscErrorCode RGGetComplement(RG,PetscBool*);
64: SLEPC_EXTERN PetscErrorCode RGSetScale(RG,PetscReal);
65: SLEPC_EXTERN PetscErrorCode RGGetScale(RG,PetscReal*);
66: SLEPC_EXTERN PetscErrorCode RGPushScale(RG,PetscReal);
67: SLEPC_EXTERN PetscErrorCode RGPopScale(RG);
68: SLEPC_EXTERN PetscErrorCode RGCheckInside(RG,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);
69: SLEPC_EXTERN PetscErrorCode RGIsAxisymmetric(RG,PetscBool,PetscBool*);
70: SLEPC_EXTERN PetscErrorCode RGCanUseConjugates(RG,PetscBool,PetscBool*);
71: SLEPC_EXTERN PetscErrorCode RGComputeContour(RG,PetscInt,PetscScalar*,PetscScalar*);
72: SLEPC_EXTERN PetscErrorCode RGComputeBoundingBox(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
73: SLEPC_EXTERN PetscErrorCode RGComputeQuadrature(RG,RGQuadRule,PetscInt,PetscScalar*,PetscScalar*,PetscScalar*);
75: SLEPC_EXTERN PetscFunctionList RGList;
76: SLEPC_EXTERN PetscErrorCode RGRegister(const char[],PetscErrorCode(*)(RG));
78: /* --------- options specific to particular regions -------- */
80: SLEPC_EXTERN PetscErrorCode RGEllipseSetParameters(RG,PetscScalar,PetscReal,PetscReal);
81: SLEPC_EXTERN PetscErrorCode RGEllipseGetParameters(RG,PetscScalar*,PetscReal*,PetscReal*);
83: SLEPC_EXTERN PetscErrorCode RGIntervalSetEndpoints(RG,PetscReal,PetscReal,PetscReal,PetscReal);
84: SLEPC_EXTERN PetscErrorCode RGIntervalGetEndpoints(RG,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
86: SLEPC_EXTERN PetscErrorCode RGPolygonSetVertices(RG,PetscInt,PetscScalar*,PetscScalar*);
87: SLEPC_EXTERN PetscErrorCode RGPolygonGetVertices(RG,PetscInt*,PetscScalar**,PetscScalar**);
89: SLEPC_EXTERN PetscErrorCode RGRingSetParameters(RG,PetscScalar,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
90: SLEPC_EXTERN PetscErrorCode RGRingGetParameters(RG,PetscScalar*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*);