Actual source code: slepcds.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 direct solver object in SLEPc
12: */
14: #if !defined(SLEPCDS_H)
15: #define SLEPCDS_H
16: #include <slepcsc.h>
17: #include <slepcfn.h>
18: #include <slepcrg.h>
20: #define DS_MAX_SOLVE 6
22: SLEPC_EXTERN PetscErrorCode DSInitializePackage(void);
23: /*S
24: DS - Direct solver (or dense system), to represent low-dimensional
25: eigenproblems that must be solved within iterative solvers. This is an
26: auxiliary object and is not normally needed by application programmers.
28: Level: beginner
30: .seealso: DSCreate()
31: S*/
32: typedef struct _p_DS* DS;
34: /*J
35: DSType - String with the name of the type of direct solver. Roughly,
36: there are as many types as problem types are available within SLEPc.
38: Level: advanced
40: .seealso: DSSetType(), DS
41: J*/
42: typedef const char* DSType;
43: #define DSHEP "hep"
44: #define DSNHEP "nhep"
45: #define DSGHEP "ghep"
46: #define DSGHIEP "ghiep"
47: #define DSGNHEP "gnhep"
48: #define DSNHEPTS "nhepts"
49: #define DSSVD "svd"
50: #define DSGSVD "gsvd"
51: #define DSPEP "pep"
52: #define DSNEP "nep"
54: /* Logging support */
55: SLEPC_EXTERN PetscClassId DS_CLASSID;
57: /*E
58: DSStateType - Indicates in which state the direct solver is
60: Level: advanced
62: .seealso: DSSetState()
63: E*/
64: typedef enum { DS_STATE_RAW,
65: DS_STATE_INTERMEDIATE,
66: DS_STATE_CONDENSED,
67: DS_STATE_TRUNCATED } DSStateType;
68: SLEPC_EXTERN const char *DSStateTypes[];
70: /*E
71: DSMatType - Used to refer to one of the matrices stored internally in DS
73: Notes:
74: The matrices preferently refer to
75: + DS_MAT_A - first matrix of eigenproblem/singular value problem
76: . DS_MAT_B - second matrix of a generalized eigenproblem
77: . DS_MAT_C - third matrix of a quadratic eigenproblem
78: . DS_MAT_T - tridiagonal matrix
79: . DS_MAT_D - diagonal matrix
80: . DS_MAT_Q - orthogonal matrix of (right) Schur vectors
81: . DS_MAT_Z - orthogonal matrix of left Schur vectors
82: . DS_MAT_X - right eigenvectors
83: . DS_MAT_Y - left eigenvectors
84: . DS_MAT_U - left singular vectors
85: . DS_MAT_V - right singular vectors
86: . DS_MAT_W - workspace matrix
87: - DS_MAT_Ex - extra matrices (x=0,..,9)
89: All matrices can have space to hold ld x ld elements, except for
90: DS_MAT_T that has space for 3 x ld elements (ld = leading dimension)
91: and DS_MAT_D that has space for just ld elements.
93: In DSPEP problems, matrices A, B, W can have space for d*ld x d*ld,
94: where d is the polynomial degree, and X can have ld x d*ld.
96: Level: advanced
98: .seealso: DSAllocate(), DSGetArray(), DSGetArrayReal(), DSVectors()
99: E*/
100: typedef enum { DS_MAT_A,
101: DS_MAT_B,
102: DS_MAT_C,
103: DS_MAT_T,
104: DS_MAT_D,
105: DS_MAT_Q,
106: DS_MAT_Z,
107: DS_MAT_X,
108: DS_MAT_Y,
109: DS_MAT_U,
110: DS_MAT_V,
111: DS_MAT_W,
112: DS_MAT_E0,
113: DS_MAT_E1,
114: DS_MAT_E2,
115: DS_MAT_E3,
116: DS_MAT_E4,
117: DS_MAT_E5,
118: DS_MAT_E6,
119: DS_MAT_E7,
120: DS_MAT_E8,
121: DS_MAT_E9,
122: DS_NUM_MAT } DSMatType;
124: /* Convenience for indexing extra matrices */
125: SLEPC_EXTERN DSMatType DSMatExtra[];
126: #define DS_NUM_EXTRA 10
128: /*E
129: DSParallelType - Indicates the parallel mode that the direct solver will use
131: Level: advanced
133: .seealso: DSSetParallel()
134: E*/
135: typedef enum { DS_PARALLEL_REDUNDANT,
136: DS_PARALLEL_SYNCHRONIZED,
137: DS_PARALLEL_DISTRIBUTED } DSParallelType;
138: SLEPC_EXTERN const char *DSParallelTypes[];
140: SLEPC_EXTERN PetscErrorCode DSCreate(MPI_Comm,DS*);
141: SLEPC_EXTERN PetscErrorCode DSSetType(DS,DSType);
142: SLEPC_EXTERN PetscErrorCode DSGetType(DS,DSType*);
143: SLEPC_EXTERN PetscErrorCode DSSetOptionsPrefix(DS,const char *);
144: SLEPC_EXTERN PetscErrorCode DSAppendOptionsPrefix(DS,const char *);
145: SLEPC_EXTERN PetscErrorCode DSGetOptionsPrefix(DS,const char *[]);
146: SLEPC_EXTERN PetscErrorCode DSSetFromOptions(DS);
147: SLEPC_EXTERN PetscErrorCode DSView(DS,PetscViewer);
148: SLEPC_EXTERN PetscErrorCode DSViewFromOptions(DS,PetscObject,const char[]);
149: SLEPC_EXTERN PetscErrorCode DSViewMat(DS,PetscViewer,DSMatType);
150: SLEPC_EXTERN PetscErrorCode DSDestroy(DS*);
151: SLEPC_EXTERN PetscErrorCode DSReset(DS);
152: SLEPC_EXTERN PetscErrorCode DSDuplicate(DS,DS*);
154: SLEPC_EXTERN PetscErrorCode DSAllocate(DS,PetscInt);
155: SLEPC_EXTERN PetscErrorCode DSGetLeadingDimension(DS,PetscInt*);
156: SLEPC_EXTERN PetscErrorCode DSSetState(DS,DSStateType);
157: SLEPC_EXTERN PetscErrorCode DSGetState(DS,DSStateType*);
158: SLEPC_EXTERN PetscErrorCode DSSetDimensions(DS,PetscInt,PetscInt,PetscInt);
159: SLEPC_EXTERN PetscErrorCode DSGetDimensions(DS,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
160: SLEPC_EXTERN PetscErrorCode DSSetBlockSize(DS,PetscInt);
161: SLEPC_EXTERN PetscErrorCode DSGetBlockSize(DS,PetscInt*);
162: SLEPC_EXTERN PetscErrorCode DSGetTruncateSize(DS,PetscInt,PetscInt,PetscInt*);
163: SLEPC_EXTERN PetscErrorCode DSTruncate(DS,PetscInt,PetscBool);
164: SLEPC_EXTERN PetscErrorCode DSSetIdentity(DS,DSMatType);
165: SLEPC_EXTERN PetscErrorCode DSSetMethod(DS,PetscInt);
166: SLEPC_EXTERN PetscErrorCode DSGetMethod(DS,PetscInt*);
167: SLEPC_EXTERN PetscErrorCode DSSetParallel(DS,DSParallelType);
168: SLEPC_EXTERN PetscErrorCode DSGetParallel(DS,DSParallelType*);
169: SLEPC_EXTERN PetscErrorCode DSSetCompact(DS,PetscBool);
170: SLEPC_EXTERN PetscErrorCode DSGetCompact(DS,PetscBool*);
171: SLEPC_EXTERN PetscErrorCode DSSetExtraRow(DS,PetscBool);
172: SLEPC_EXTERN PetscErrorCode DSGetExtraRow(DS,PetscBool*);
173: SLEPC_EXTERN PetscErrorCode DSSetRefined(DS,PetscBool);
174: SLEPC_EXTERN PetscErrorCode DSGetRefined(DS,PetscBool*);
175: SLEPC_EXTERN PetscErrorCode DSGetMat(DS,DSMatType,Mat*);
176: SLEPC_EXTERN PetscErrorCode DSRestoreMat(DS,DSMatType,Mat*);
177: SLEPC_EXTERN PetscErrorCode DSGetArray(DS,DSMatType,PetscScalar*[]);
178: SLEPC_EXTERN PetscErrorCode DSRestoreArray(DS,DSMatType,PetscScalar*[]);
179: SLEPC_EXTERN PetscErrorCode DSGetArrayReal(DS,DSMatType,PetscReal*[]);
180: SLEPC_EXTERN PetscErrorCode DSRestoreArrayReal(DS,DSMatType,PetscReal*[]);
181: SLEPC_EXTERN PetscErrorCode DSVectors(DS,DSMatType,PetscInt*,PetscReal*);
182: SLEPC_EXTERN PetscErrorCode DSSolve(DS,PetscScalar*,PetscScalar*);
183: SLEPC_EXTERN PetscErrorCode DSSort(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
184: SLEPC_EXTERN PetscErrorCode DSSortWithPermutation(DS,PetscInt*,PetscScalar*,PetscScalar*);
185: SLEPC_EXTERN PetscErrorCode DSSynchronize(DS,PetscScalar*,PetscScalar*);
186: SLEPC_EXTERN PetscErrorCode DSCopyMat(DS,DSMatType,PetscInt,PetscInt,Mat,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
187: SLEPC_EXTERN PetscErrorCode DSMatGetSize(DS,DSMatType,PetscInt*,PetscInt*);
188: SLEPC_EXTERN PetscErrorCode DSMatIsHermitian(DS,DSMatType,PetscBool*);
189: SLEPC_EXTERN PetscErrorCode DSSetSlepcSC(DS,SlepcSC);
190: SLEPC_EXTERN PetscErrorCode DSGetSlepcSC(DS,SlepcSC*);
191: SLEPC_EXTERN PetscErrorCode DSUpdateExtraRow(DS);
192: SLEPC_EXTERN PetscErrorCode DSCond(DS,PetscReal*);
193: SLEPC_EXTERN PetscErrorCode DSTranslateHarmonic(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
194: SLEPC_EXTERN PetscErrorCode DSTranslateRKS(DS,PetscScalar);
195: SLEPC_EXTERN PetscErrorCode DSOrthogonalize(DS,DSMatType,PetscInt,PetscInt*);
196: SLEPC_EXTERN PetscErrorCode DSPseudoOrthogonalize(DS,DSMatType,PetscInt,PetscReal*,PetscInt*,PetscReal*);
198: /* --------- options specific to particular solvers -------- */
200: SLEPC_EXTERN PetscErrorCode DSSVDSetDimensions(DS,PetscInt);
201: SLEPC_EXTERN PetscErrorCode DSSVDGetDimensions(DS,PetscInt*);
202: SLEPC_EXTERN PetscErrorCode DSGSVDSetDimensions(DS,PetscInt,PetscInt);
203: SLEPC_EXTERN PetscErrorCode DSGSVDGetDimensions(DS,PetscInt*,PetscInt*);
205: SLEPC_EXTERN PetscErrorCode DSPEPSetDegree(DS,PetscInt);
206: SLEPC_EXTERN PetscErrorCode DSPEPGetDegree(DS,PetscInt*);
207: SLEPC_EXTERN PetscErrorCode DSPEPSetCoefficients(DS,PetscReal*);
208: SLEPC_EXTERN PetscErrorCode DSPEPGetCoefficients(DS,PetscReal**);
210: SLEPC_EXTERN PetscErrorCode DSNEPSetFN(DS,PetscInt,FN*);
211: SLEPC_EXTERN PetscErrorCode DSNEPGetFN(DS,PetscInt,FN*);
212: SLEPC_EXTERN PetscErrorCode DSNEPGetNumFN(DS,PetscInt*);
213: SLEPC_EXTERN PetscErrorCode DSNEPSetMinimality(DS,PetscInt);
214: SLEPC_EXTERN PetscErrorCode DSNEPGetMinimality(DS,PetscInt*);
215: SLEPC_EXTERN PetscErrorCode DSNEPSetRefine(DS,PetscReal,PetscInt);
216: SLEPC_EXTERN PetscErrorCode DSNEPGetRefine(DS,PetscReal*,PetscInt*);
217: SLEPC_EXTERN PetscErrorCode DSNEPSetIntegrationPoints(DS,PetscInt);
218: SLEPC_EXTERN PetscErrorCode DSNEPGetIntegrationPoints(DS,PetscInt*);
219: SLEPC_EXTERN PetscErrorCode DSNEPSetSamplingSize(DS,PetscInt);
220: SLEPC_EXTERN PetscErrorCode DSNEPGetSamplingSize(DS,PetscInt*);
221: SLEPC_EXTERN PetscErrorCode DSNEPSetRG(DS,RG);
222: SLEPC_EXTERN PetscErrorCode DSNEPGetRG(DS,RG*);
223: SLEPC_EXTERN PetscErrorCode DSNEPSetComputeMatrixFunction(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*);
224: SLEPC_EXTERN PetscErrorCode DSNEPGetComputeMatrixFunction(DS,PetscErrorCode (**)(DS,PetscScalar,PetscBool,DSMatType,void*),void**);
226: SLEPC_EXTERN PetscFunctionList DSList;
227: SLEPC_EXTERN PetscErrorCode DSRegister(const char[],PetscErrorCode(*)(DS));
229: #endif