Actual source code: slepcsvd.h
slepc-3.19.2 2023-09-05
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 SLEPc's singular value solvers
12: */
14: #if !defined(SLEPCSVD_H)
15: #define SLEPCSVD_H
17: #include <slepceps.h>
18: #include <slepcbv.h>
19: #include <slepcds.h>
21: /* SUBMANSEC = SVD */
23: SLEPC_EXTERN PetscErrorCode SVDInitializePackage(void);
25: /*S
26: SVD - Abstract SLEPc object that manages all the singular value
27: problem solvers.
29: Level: beginner
31: .seealso: SVDCreate()
32: S*/
33: typedef struct _p_SVD* SVD;
35: /*J
36: SVDType - String with the name of a SLEPc singular value solver
38: Level: beginner
40: .seealso: SVDSetType(), SVD
41: J*/
42: typedef const char* SVDType;
43: #define SVDCROSS "cross"
44: #define SVDCYCLIC "cyclic"
45: #define SVDLAPACK "lapack"
46: #define SVDLANCZOS "lanczos"
47: #define SVDTRLANCZOS "trlanczos"
48: #define SVDRANDOMIZED "randomized"
49: #define SVDSCALAPACK "scalapack"
50: #define SVDKSVD "ksvd"
51: #define SVDELEMENTAL "elemental"
52: #define SVDPRIMME "primme"
54: /* Logging support */
55: SLEPC_EXTERN PetscClassId SVD_CLASSID;
57: /*E
58: SVDProblemType - Determines the type of singular value problem
60: Level: beginner
62: .seealso: SVDSetProblemType(), SVDGetProblemType()
63: E*/
64: typedef enum { SVD_STANDARD=1,
65: SVD_GENERALIZED, /* GSVD */
66: SVD_HYPERBOLIC /* HSVD */
67: } SVDProblemType;
69: /*E
70: SVDWhich - Determines whether largest or smallest singular triplets
71: are to be computed
73: Level: intermediate
75: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
76: E*/
77: typedef enum { SVD_LARGEST,
78: SVD_SMALLEST } SVDWhich;
80: /*E
81: SVDErrorType - The error type used to assess accuracy of computed solutions
83: Level: intermediate
85: .seealso: SVDComputeError()
86: E*/
87: typedef enum { SVD_ERROR_ABSOLUTE,
88: SVD_ERROR_RELATIVE,
89: SVD_ERROR_NORM } SVDErrorType;
90: SLEPC_EXTERN const char *SVDErrorTypes[];
92: /*E
93: SVDConv - Determines the convergence test
95: Level: intermediate
97: .seealso: SVDSetConvergenceTest(), SVDSetConvergenceTestFunction()
98: E*/
99: typedef enum { SVD_CONV_ABS,
100: SVD_CONV_REL,
101: SVD_CONV_NORM,
102: SVD_CONV_MAXIT,
103: SVD_CONV_USER } SVDConv;
105: /*E
106: SVDStop - Determines the stopping test
108: Level: advanced
110: .seealso: SVDSetStoppingTest(), SVDSetStoppingTestFunction()
111: E*/
112: typedef enum { SVD_STOP_BASIC,
113: SVD_STOP_USER } SVDStop;
115: /*E
116: SVDConvergedReason - Reason a singular value solver was said to
117: have converged or diverged
119: Level: intermediate
121: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
122: E*/
123: typedef enum {/* converged */
124: SVD_CONVERGED_TOL = 1,
125: SVD_CONVERGED_USER = 2,
126: SVD_CONVERGED_MAXIT = 3,
127: /* diverged */
128: SVD_DIVERGED_ITS = -1,
129: SVD_DIVERGED_BREAKDOWN = -2,
130: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
131: SLEPC_EXTERN const char *const*SVDConvergedReasons;
133: SLEPC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
134: SLEPC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
135: SLEPC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
136: SLEPC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
137: SLEPC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
138: SLEPC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
139: SLEPC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
140: SLEPC_EXTERN PetscErrorCode SVDSetProblemType(SVD,SVDProblemType);
141: SLEPC_EXTERN PetscErrorCode SVDGetProblemType(SVD,SVDProblemType*);
142: SLEPC_EXTERN PetscErrorCode SVDIsGeneralized(SVD,PetscBool*);
143: SLEPC_EXTERN PetscErrorCode SVDIsHyperbolic(SVD,PetscBool*);
144: SLEPC_EXTERN PetscErrorCode SVDSetOperators(SVD,Mat,Mat);
145: PETSC_DEPRECATED_FUNCTION("Use SVDSetOperators()") static inline PetscErrorCode SVDSetOperator(SVD svd,Mat A) {return SVDSetOperators(svd,A,PETSC_NULLPTR);}
146: SLEPC_EXTERN PetscErrorCode SVDGetOperators(SVD,Mat*,Mat*);
147: PETSC_DEPRECATED_FUNCTION("Use SVDGetOperators()") static inline PetscErrorCode SVDGetOperator(SVD svd,Mat *A) {return SVDGetOperators(svd,A,PETSC_NULLPTR);}
148: SLEPC_EXTERN PetscErrorCode SVDSetSignature(SVD,Vec);
149: SLEPC_EXTERN PetscErrorCode SVDGetSignature(SVD,Vec*);
150: SLEPC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec[],PetscInt,Vec[]);
151: PETSC_DEPRECATED_FUNCTION("Use SVDSetInitialSpaces()") static inline PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt nr,Vec *isr) {return SVDSetInitialSpaces(svd,nr,isr,0,PETSC_NULLPTR);}
152: PETSC_DEPRECATED_FUNCTION("Use SVDSetInitialSpaces()") static inline PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt nl,Vec *isl) {return SVDSetInitialSpaces(svd,0,PETSC_NULLPTR,nl,isl);}
153: SLEPC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
154: SLEPC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
155: SLEPC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
156: SLEPC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
157: SLEPC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
158: SLEPC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
159: SLEPC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
160: SLEPC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
161: SLEPC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
162: SLEPC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
163: SLEPC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
164: SLEPC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
165: SLEPC_EXTERN PetscErrorCode SVDSetDSType(SVD);
166: SLEPC_EXTERN PetscErrorCode SVDSetUp(SVD);
167: SLEPC_EXTERN PetscErrorCode SVDSolve(SVD);
168: SLEPC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
169: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,PetscErrorCode (*)(SVD,PetscReal,PetscReal,PetscReal*,void*),void*,PetscErrorCode (*)(void*));
170: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
171: SLEPC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
172: SLEPC_EXTERN PetscErrorCode SVDConvergedAbsolute(SVD,PetscReal,PetscReal,PetscReal*,void*);
173: SLEPC_EXTERN PetscErrorCode SVDConvergedRelative(SVD,PetscReal,PetscReal,PetscReal*,void*);
174: SLEPC_EXTERN PetscErrorCode SVDConvergedNorm(SVD,PetscReal,PetscReal,PetscReal*,void*);
175: SLEPC_EXTERN PetscErrorCode SVDConvergedMaxIt(SVD,PetscReal,PetscReal,PetscReal*,void*);
176: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
177: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
178: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
179: SLEPC_EXTERN PetscErrorCode SVDStoppingBasic(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
180: SLEPC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
181: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
182: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
183: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
184: PETSC_DEPRECATED_FUNCTION("Use SVDComputeError()") static inline PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *r) {return SVDComputeError(svd,i,SVD_ERROR_RELATIVE,r);}
185: PETSC_DEPRECATED_FUNCTION("Use SVDComputeError() with SVD_ERROR_ABSOLUTE") static inline PetscErrorCode SVDComputeResidualNorms(SVD svd,PetscInt i,PetscReal *r1,PETSC_UNUSED PetscReal *r2) {return SVDComputeError(svd,i,SVD_ERROR_ABSOLUTE,r1);}
186: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
187: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
188: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
189: PETSC_DEPRECATED_FUNCTION("Use SVDErrorView()") static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
190: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
191: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
192: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
193: PETSC_DEPRECATED_FUNCTION("Use SVDConvergedReasonView() (since version 3.14)") static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
194: PETSC_DEPRECATED_FUNCTION("Use SVDConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
195: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
196: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
197: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
198: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
199: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
200: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
201: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
202: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
203: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
205: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
206: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
207: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
208: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);
210: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
211: SLEPC_EXTERN PetscErrorCode SVDMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
212: SLEPC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
213: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
214: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
215: SLEPC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
216: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
217: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
218: SLEPC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
219: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
220: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
221: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
222: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat**);
223: SLEPC_EXTERN PetscErrorCode SVDMonitorConditioning(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
225: SLEPC_EXTERN PetscFunctionList SVDList;
226: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
227: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
228: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
229: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
230: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));
232: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
234: /* --------- options specific to particular solvers -------- */
236: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
237: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
238: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
239: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
241: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
242: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
243: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
244: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
246: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
247: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
249: /*E
250: SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
251: TRLanczos GSVD solver
253: Level: advanced
255: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
256: E*/
257: typedef enum {
258: SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
259: SVD_TRLANCZOS_GBIDIAG_UPPER, /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
260: SVD_TRLANCZOS_GBIDIAG_LOWER /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
261: } SVDTRLanczosGBidiag;
262: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];
264: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
265: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
266: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
267: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
268: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
269: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
270: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
271: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
272: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
273: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
274: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
275: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
276: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
277: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);
279: /*E
280: SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library
282: Level: advanced
284: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
285: E*/
286: typedef enum { SVD_PRIMME_HYBRID=1,
287: SVD_PRIMME_NORMALEQUATIONS,
288: SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
289: SLEPC_EXTERN const char *SVDPRIMMEMethods[];
291: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
292: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
293: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
294: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);
296: /*E
297: SVDKSVDEigenMethod - determines the method to solve the eigenproblem within KSVD
299: Level: advanced
301: .seealso: SVDKSVDSetEigenMethod(), SVDKSVDGetEigenMethod()
302: E*/
303: typedef enum { SVD_KSVD_EIGEN_MRRR=1,
304: SVD_KSVD_EIGEN_DC,
305: SVD_KSVD_EIGEN_ELPA } SVDKSVDEigenMethod;
306: SLEPC_EXTERN const char *SVDKSVDEigenMethods[];
308: /*E
309: SVDKSVDPolarMethod - determines the method to compute the polar decomposition within KSVD
311: Level: advanced
313: .seealso: SVDKSVDSetPolarMethod(), SVDKSVDGetPolarMethod()
314: E*/
315: typedef enum { SVD_KSVD_POLAR_QDWH=1,
316: SVD_KSVD_POLAR_ZOLOPD } SVDKSVDPolarMethod;
317: SLEPC_EXTERN const char *SVDKSVDPolarMethods[];
319: SLEPC_EXTERN PetscErrorCode SVDKSVDSetEigenMethod(SVD,SVDKSVDEigenMethod);
320: SLEPC_EXTERN PetscErrorCode SVDKSVDGetEigenMethod(SVD,SVDKSVDEigenMethod*);
321: SLEPC_EXTERN PetscErrorCode SVDKSVDSetPolarMethod(SVD,SVDKSVDPolarMethod);
322: SLEPC_EXTERN PetscErrorCode SVDKSVDGetPolarMethod(SVD,SVDKSVDPolarMethod*);
324: #endif