Actual source code: slepcsvd.h
slepc-main 2025-03-25
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: #pragma once
16: #include <slepceps.h>
17: #include <slepcbv.h>
18: #include <slepcds.h>
20: /* SUBMANSEC = SVD */
22: SLEPC_EXTERN PetscErrorCode SVDInitializePackage(void);
23: SLEPC_EXTERN PetscErrorCode SVDFinalizePackage(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 = 2, /* GSVD */
66: SVD_HYPERBOLIC = 3 /* 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,
114: SVD_STOP_THRESHOLD } SVDStop;
116: /*E
117: SVDConvergedReason - Reason a singular value solver was said to
118: have converged or diverged
120: Level: intermediate
122: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
123: E*/
124: typedef enum {/* converged */
125: SVD_CONVERGED_TOL = 1,
126: SVD_CONVERGED_USER = 2,
127: SVD_CONVERGED_MAXIT = 3,
128: /* diverged */
129: SVD_DIVERGED_ITS = -1,
130: SVD_DIVERGED_BREAKDOWN = -2,
131: SVD_DIVERGED_SYMMETRY_LOST = -3,
132: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
133: SLEPC_EXTERN const char *const*SVDConvergedReasons;
135: /*S
136: SVDStoppingCtx - Data structure (C struct) to hold additional information to
137: be used in some stopping test functions.
139: Level: advanced
141: .seealso: SVDSetStoppingTestFunction()
142: S*/
143: struct _n_SVDStoppingCtx {
144: PetscReal firstsv; /* the value of the first converged singular value */
145: PetscReal lastsv; /* the value of the last converged singular value */
146: PetscReal thres; /* threshold set with SVDSetThreshold() */
147: PetscBool threlative; /* threshold is relative */
148: SVDWhich which; /* which singular values are being computed */
149: };
150: typedef struct _n_SVDStoppingCtx* SVDStoppingCtx;
152: SLEPC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
153: SLEPC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
154: SLEPC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
155: SLEPC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
156: SLEPC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
157: SLEPC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
158: SLEPC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
159: SLEPC_EXTERN PetscErrorCode SVDSetProblemType(SVD,SVDProblemType);
160: SLEPC_EXTERN PetscErrorCode SVDGetProblemType(SVD,SVDProblemType*);
161: SLEPC_EXTERN PetscErrorCode SVDIsGeneralized(SVD,PetscBool*);
162: SLEPC_EXTERN PetscErrorCode SVDIsHyperbolic(SVD,PetscBool*);
163: SLEPC_EXTERN PetscErrorCode SVDSetOperators(SVD,Mat,Mat);
164: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDSetOperators()", ) static inline PetscErrorCode SVDSetOperator(SVD svd,Mat A) {return SVDSetOperators(svd,A,PETSC_NULLPTR);}
165: SLEPC_EXTERN PetscErrorCode SVDGetOperators(SVD,Mat*,Mat*);
166: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDGetOperators()", ) static inline PetscErrorCode SVDGetOperator(SVD svd,Mat *A) {return SVDGetOperators(svd,A,PETSC_NULLPTR);}
167: SLEPC_EXTERN PetscErrorCode SVDSetSignature(SVD,Vec);
168: SLEPC_EXTERN PetscErrorCode SVDGetSignature(SVD,Vec);
169: SLEPC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec[],PetscInt,Vec[]);
170: PETSC_DEPRECATED_FUNCTION(3, 1, 0, "SVDSetInitialSpaces()", ) static inline PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt nr,Vec *isr) {return SVDSetInitialSpaces(svd,nr,isr,0,PETSC_NULLPTR);}
171: PETSC_DEPRECATED_FUNCTION(3, 1, 0, "SVDSetInitialSpaces()", ) static inline PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt nl,Vec *isl) {return SVDSetInitialSpaces(svd,0,PETSC_NULLPTR,nl,isl);}
172: SLEPC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
173: SLEPC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
174: SLEPC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
175: SLEPC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
176: SLEPC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
177: SLEPC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
178: SLEPC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
179: SLEPC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
180: SLEPC_EXTERN PetscErrorCode SVDSetThreshold(SVD,PetscReal,PetscBool);
181: SLEPC_EXTERN PetscErrorCode SVDGetThreshold(SVD,PetscReal*,PetscBool*);
182: SLEPC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
183: SLEPC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
184: SLEPC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
185: SLEPC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
186: SLEPC_EXTERN PetscErrorCode SVDSetDSType(SVD);
187: SLEPC_EXTERN PetscErrorCode SVDSetUp(SVD);
188: SLEPC_EXTERN PetscErrorCode SVDSolve(SVD);
189: SLEPC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
190: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
191: SLEPC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
192: SLEPC_EXTERN PetscErrorCode SVDConvergedAbsolute(SVD,PetscReal,PetscReal,PetscReal*,void*);
193: SLEPC_EXTERN PetscErrorCode SVDConvergedRelative(SVD,PetscReal,PetscReal,PetscReal*,void*);
194: SLEPC_EXTERN PetscErrorCode SVDConvergedNorm(SVD,PetscReal,PetscReal,PetscReal*,void*);
195: SLEPC_EXTERN PetscErrorCode SVDConvergedMaxIt(SVD,PetscReal,PetscReal,PetscReal*,void*);
196: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
197: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
198: SLEPC_EXTERN PetscErrorCode SVDStoppingBasic(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
199: SLEPC_EXTERN PetscErrorCode SVDStoppingThreshold(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
200: SLEPC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
201: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
202: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
203: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
204: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "SVDComputeError()", ) static inline PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *r) {return SVDComputeError(svd,i,SVD_ERROR_RELATIVE,r);}
205: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "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);}
206: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
207: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
208: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
209: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "SVDErrorView()", ) static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
210: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
211: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
212: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
213: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonView()", ) static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
214: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
215: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
216: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
217: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
218: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
219: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
220: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
221: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
222: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
223: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
225: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
226: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscCtxDestroyFn*);
227: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
228: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);
230: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
231: SLEPC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
232: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
233: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
234: SLEPC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
235: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
236: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
237: SLEPC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
238: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
239: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
240: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
241: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat**);
242: SLEPC_EXTERN PetscErrorCode SVDMonitorConditioning(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
244: SLEPC_EXTERN PetscFunctionList SVDList;
245: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
246: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
247: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
248: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
249: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));
251: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
252: SLEPC_EXTERN PetscErrorCode SVDReallocateSolution(SVD,PetscInt);
254: /*S
255: SVDConvergenceTestFn - A prototype of an SVD convergence test function that would be passed to SVDSetConvergenceTestFunction()
257: Calling Sequence:
258: + svd - singular value solver context obtained from SVDCreate()
259: . sigma - computed singular value
260: . res - residual norm associated to the singular triplet
261: . errest - [output] computed error estimate
262: - ctx - [optional] user-defined context for private data for the
263: convergence test routine (may be NULL)
265: Level: advanced
267: .seealso: SVDSetConvergenceTestFunction()
268: S*/
269: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SVDConvergenceTestFn)(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx);
271: /*S
272: SVDStoppingTestFn - A prototype of an SVD stopping test function that would be passed to SVDSetStoppingTestFunction()
274: Calling Sequence:
275: + svd - singular value solver context obtained from SVDCreate()
276: . its - current number of iterations
277: . max_it - maximum number of iterations
278: . nconv - number of currently converged singular triplets
279: . nsv - number of requested singular triplets
280: . reason - [output] result of the stopping test
281: - ctx - [optional] user-defined context for private data for the
282: stopping test routine (may be NULL)
284: Level: advanced
286: .seealso: SVDSetStoppingTestFunction()
287: S*/
288: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SVDStoppingTestFn)(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx);
290: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,SVDConvergenceTestFn*,void*,PetscCtxDestroyFn*);
291: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,SVDStoppingTestFn*,void*,PetscCtxDestroyFn*);
293: /* --------- options specific to particular solvers -------- */
295: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
296: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
297: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
298: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
300: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
301: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
302: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
303: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
305: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
306: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
308: /*E
309: SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
310: TRLanczos GSVD solver
312: Level: advanced
314: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
315: E*/
316: typedef enum {
317: SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
318: SVD_TRLANCZOS_GBIDIAG_UPPER, /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
319: SVD_TRLANCZOS_GBIDIAG_LOWER /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
320: } SVDTRLanczosGBidiag;
321: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];
323: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
324: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
325: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
326: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
327: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
328: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
329: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
330: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
331: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
332: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
333: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
334: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
335: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
336: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);
338: /*E
339: SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library
341: Level: advanced
343: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
344: E*/
345: typedef enum { SVD_PRIMME_HYBRID = 1,
346: SVD_PRIMME_NORMALEQUATIONS = 2,
347: SVD_PRIMME_AUGMENTED = 3 } SVDPRIMMEMethod;
348: SLEPC_EXTERN const char *SVDPRIMMEMethods[];
350: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
351: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
352: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
353: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);
355: /*E
356: SVDKSVDEigenMethod - determines the method to solve the eigenproblem within KSVD
358: Level: advanced
360: .seealso: SVDKSVDSetEigenMethod(), SVDKSVDGetEigenMethod()
361: E*/
362: typedef enum { SVD_KSVD_EIGEN_MRRR = 1,
363: SVD_KSVD_EIGEN_DC = 2,
364: SVD_KSVD_EIGEN_ELPA = 3 } SVDKSVDEigenMethod;
365: SLEPC_EXTERN const char *SVDKSVDEigenMethods[];
367: /*E
368: SVDKSVDPolarMethod - determines the method to compute the polar decomposition within KSVD
370: Level: advanced
372: .seealso: SVDKSVDSetPolarMethod(), SVDKSVDGetPolarMethod()
373: E*/
374: typedef enum { SVD_KSVD_POLAR_QDWH = 1,
375: SVD_KSVD_POLAR_ZOLOPD = 2 } SVDKSVDPolarMethod;
376: SLEPC_EXTERN const char *SVDKSVDPolarMethods[];
378: SLEPC_EXTERN PetscErrorCode SVDKSVDSetEigenMethod(SVD,SVDKSVDEigenMethod);
379: SLEPC_EXTERN PetscErrorCode SVDKSVDGetEigenMethod(SVD,SVDKSVDEigenMethod*);
380: SLEPC_EXTERN PetscErrorCode SVDKSVDSetPolarMethod(SVD,SVDKSVDPolarMethod);
381: SLEPC_EXTERN PetscErrorCode SVDKSVDGetPolarMethod(SVD,SVDKSVDPolarMethod*);