Actual source code: slepcsvd.h
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 - SLEPc object that manages all the singular value problem solvers.
28: Level: beginner
30: .seealso: [](ch:svd), `SVDCreate()`
31: S*/
32: typedef struct _p_SVD* SVD;
34: /*J
35: SVDType - String with the name of a singular value solver.
37: Level: beginner
39: .seealso: [](ch:svd), `SVDSetType()`, `SVD`
40: J*/
41: typedef const char *SVDType;
42: #define SVDCROSS "cross"
43: #define SVDCYCLIC "cyclic"
44: #define SVDLANCZOS "lanczos"
45: #define SVDTRLANCZOS "trlanczos"
46: #define SVDRANDOMIZED "randomized"
47: #define SVDLAPACK "lapack"
48: #define SVDSCALAPACK "scalapack"
49: #define SVDKSVD "ksvd"
50: #define SVDELEMENTAL "elemental"
51: #define SVDPRIMME "primme"
53: /* Logging support */
54: SLEPC_EXTERN PetscClassId SVD_CLASSID;
56: /*E
57: SVDProblemType - Determines the type of singular value problem
59: Level: beginner
61: .seealso: [](ch:svd), `SVDSetProblemType()`, `SVDGetProblemType()`
62: E*/
63: typedef enum { SVD_STANDARD = 1,
64: SVD_GENERALIZED = 2, /* GSVD */
65: SVD_HYPERBOLIC = 3 /* HSVD */
66: } SVDProblemType;
68: /*E
69: SVDWhich - Determines whether largest or smallest singular triplets
70: are to be computed
72: Level: intermediate
74: .seealso: [](ch:svd), `SVDSetWhichSingularTriplets()`, `SVDGetWhichSingularTriplets()`
75: E*/
76: typedef enum { SVD_LARGEST,
77: SVD_SMALLEST } SVDWhich;
79: /*E
80: SVDErrorType - The error type used to assess accuracy of computed solutions
82: Level: intermediate
84: .seealso: [](ch:svd), `SVDComputeError()`
85: E*/
86: typedef enum { SVD_ERROR_ABSOLUTE,
87: SVD_ERROR_RELATIVE,
88: SVD_ERROR_NORM } SVDErrorType;
89: SLEPC_EXTERN const char *SVDErrorTypes[];
91: /*E
92: SVDConv - Determines the convergence test
94: Level: intermediate
96: .seealso: [](ch:svd), `SVDSetConvergenceTest()`, `SVDSetConvergenceTestFunction()`
97: E*/
98: typedef enum { SVD_CONV_ABS,
99: SVD_CONV_REL,
100: SVD_CONV_NORM,
101: SVD_CONV_MAXIT,
102: SVD_CONV_USER } SVDConv;
104: /*E
105: SVDStop - Determines the stopping test
107: Level: advanced
109: .seealso: [](ch:svd), `SVDSetStoppingTest()`, `SVDSetStoppingTestFunction()`
110: E*/
111: typedef enum { SVD_STOP_BASIC,
112: SVD_STOP_USER,
113: SVD_STOP_THRESHOLD } SVDStop;
115: /*E
116: SVDConvergedReason - Reason a singular value solver was determined to have
117: converged or diverged
119: Values:
120: + `SVD_CONVERGED_TOL` - converged up to tolerance
121: . `SVD_CONVERGED_USER` - converged due to a user-defined condition
122: . `SVD_CONVERGED_MAXIT` - reached maximum number of iterations with `SVD_CONV_MAXIT` criterion
123: . `SVD_DIVERGED_ITS` - exceeded the maximum number of allowed iterations
124: . `SVD_DIVERGED_BREAKDOWN` - generic breakdown in method
125: - `SVD_DIVERGED_SYMMETRY_LOST` - underlying indefinite eigensolver was not able to keep symmetry
127: Level: intermediate
129: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDSetTolerances()`
130: E*/
131: typedef enum {/* converged */
132: SVD_CONVERGED_TOL = 1,
133: SVD_CONVERGED_USER = 2,
134: SVD_CONVERGED_MAXIT = 3,
135: /* diverged */
136: SVD_DIVERGED_ITS = -1,
137: SVD_DIVERGED_BREAKDOWN = -2,
138: SVD_DIVERGED_SYMMETRY_LOST = -3,
139: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
140: SLEPC_EXTERN const char *const*SVDConvergedReasons;
142: /*S
143: SVDStoppingCtx - Data structure (C struct) to hold additional information to
144: be used in some stopping test functions.
146: Level: advanced
148: .seealso: [](ch:svd), `SVDSetStoppingTestFunction()`
149: S*/
150: struct _n_SVDStoppingCtx {
151: PetscReal firstsv; /* the value of the first converged singular value */
152: PetscReal lastsv; /* the value of the last converged singular value */
153: PetscReal thres; /* threshold set with SVDSetThreshold() */
154: PetscBool threlative; /* threshold is relative */
155: SVDWhich which; /* which singular values are being computed */
156: };
157: typedef struct _n_SVDStoppingCtx* SVDStoppingCtx;
159: SLEPC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
160: SLEPC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
161: SLEPC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
162: SLEPC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
163: SLEPC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
164: SLEPC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
165: SLEPC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
166: SLEPC_EXTERN PetscErrorCode SVDSetProblemType(SVD,SVDProblemType);
167: SLEPC_EXTERN PetscErrorCode SVDGetProblemType(SVD,SVDProblemType*);
168: SLEPC_EXTERN PetscErrorCode SVDIsGeneralized(SVD,PetscBool*);
169: SLEPC_EXTERN PetscErrorCode SVDIsHyperbolic(SVD,PetscBool*);
170: SLEPC_EXTERN PetscErrorCode SVDSetOperators(SVD,Mat,Mat);
171: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDSetOperators()", ) static inline PetscErrorCode SVDSetOperator(SVD svd,Mat A) {return SVDSetOperators(svd,A,PETSC_NULLPTR);}
172: SLEPC_EXTERN PetscErrorCode SVDGetOperators(SVD,Mat*,Mat*);
173: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDGetOperators()", ) static inline PetscErrorCode SVDGetOperator(SVD svd,Mat *A) {return SVDGetOperators(svd,A,PETSC_NULLPTR);}
174: SLEPC_EXTERN PetscErrorCode SVDSetSignature(SVD,Vec);
175: SLEPC_EXTERN PetscErrorCode SVDGetSignature(SVD,Vec);
176: SLEPC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec[],PetscInt,Vec[]);
177: 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);}
178: 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);}
179: SLEPC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
180: SLEPC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
181: SLEPC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
182: SLEPC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
183: SLEPC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
184: SLEPC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
185: SLEPC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
186: SLEPC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
187: SLEPC_EXTERN PetscErrorCode SVDSetThreshold(SVD,PetscReal,PetscBool);
188: SLEPC_EXTERN PetscErrorCode SVDGetThreshold(SVD,PetscReal*,PetscBool*);
189: SLEPC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
190: SLEPC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char[]);
191: SLEPC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char[]);
192: SLEPC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
193: SLEPC_EXTERN PetscErrorCode SVDSetDSType(SVD);
194: SLEPC_EXTERN PetscErrorCode SVDSetUp(SVD);
195: SLEPC_EXTERN PetscErrorCode SVDSolve(SVD);
196: SLEPC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
197: SLEPC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
198: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
199: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
200: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
201: 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);}
202: 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);}
203: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
204: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
205: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
206: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "SVDErrorView()", ) static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
207: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
208: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
209: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
210: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonView()", ) static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
211: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
212: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
213: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
214: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
215: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
216: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
217: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
218: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
219: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
220: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
222: /*S
223: SVDMonitorFn - A function prototype for functions provided to `SVDMonitorSet()`.
225: Calling Sequence:
226: + svd - the singular value solver context
227: . its - iteration number
228: . nconv - number of converged singular triplets
229: . sigma - singular values
230: . errest - relative error estimates for each singular triplet
231: . nest - number of error estimates
232: - ctx - optional monitoring context, as provided with `SVDMonitorSet()`
234: Level: intermediate
236: .seealso: [](ch:svd), `SVDMonitorSet()`
237: S*/
238: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorFn(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,void *ctx);
240: /*S
241: SVDMonitorRegisterFn - A function prototype for functions provided to `SVDMonitorRegister()`.
243: Calling Sequence:
244: + svd - the singular value solver context
245: . its - iteration number
246: . nconv - number of converged singular triplets
247: . sigma - singular values
248: . errest - relative error estimates for each singular triplet
249: . nest - number of error estimates
250: - ctx - `PetscViewerAndFormat` object
252: Level: advanced
254: Note:
255: This is an `SVDMonitorFn` specialized for a context of `PetscViewerAndFormat`.
257: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorRegister()`, `SVDMonitorFn`, `SVDMonitorRegisterCreateFn`, `SVDMonitorRegisterDestroyFn`
258: S*/
259: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterFn(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *ctx);
261: /*S
262: SVDMonitorRegisterCreateFn - A function prototype for functions that do
263: the creation when provided to `SVDMonitorRegister()`.
265: Calling Sequence:
266: + viewer - the viewer to be used with the `SVDMonitorRegisterFn`
267: . format - the format of the viewer
268: . ctx - a context for the monitor
269: - result - a `PetscViewerAndFormat` object
271: Level: advanced
273: .seealso: [](ch:svd), `SVDMonitorRegisterFn`, `SVDMonitorSet()`, `SVDMonitorRegister()`, `SVDMonitorFn`, `SVDMonitorRegisterDestroyFn`
274: S*/
275: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **result);
277: /*S
278: SVDMonitorRegisterDestroyFn - A function prototype for functions that do the after
279: use destruction when provided to `SVDMonitorRegister()`.
281: Calling Sequence:
282: . vf - a `PetscViewerAndFormat` object to be destroyed, including any context
284: Level: advanced
286: .seealso: [](ch:svd), `SVDMonitorRegisterFn`, `SVDMonitorSet()`, `SVDMonitorRegister()`, `SVDMonitorFn`, `SVDMonitorRegisterCreateFn`
287: S*/
288: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterDestroyFn(PetscViewerAndFormat **result);
290: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal[],PetscReal[],PetscInt);
291: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,SVDMonitorFn,void*,PetscCtxDestroyFn*);
292: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
293: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);
295: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
296: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorFirst;
297: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorFirstDrawLG;
298: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorFirstDrawLGCreate;
299: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorAll;
300: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorAllDrawLG;
301: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorAllDrawLGCreate;
302: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorConverged;
303: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorConvergedCreate;
304: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorConvergedDrawLG;
305: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorConvergedDrawLGCreate;
306: SLEPC_EXTERN SVDMonitorRegisterDestroyFn SVDMonitorConvergedDestroy;
307: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorConditioning;
309: SLEPC_EXTERN PetscFunctionList SVDList;
310: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
311: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
312: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
313: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
314: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,SVDMonitorRegisterFn*,SVDMonitorRegisterCreateFn*,SVDMonitorRegisterDestroyFn*);
316: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
317: SLEPC_EXTERN PetscErrorCode SVDReallocateSolution(SVD,PetscInt);
319: /*S
320: SVDConvergenceTestFn - A prototype of an SVD convergence test function that would be passed to SVDSetConvergenceTestFunction()
322: Calling Sequence:
323: + svd - the singular value solver context
324: . sigma - computed singular value
325: . res - residual norm associated to the singular triplet
326: . errest - [output] computed error estimate
327: - ctx - [optional] user-defined context for private data for the
328: convergence test routine (may be NULL)
330: Level: advanced
332: .seealso: [](ch:svd), `SVDSetConvergenceTestFunction()`
333: S*/
334: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDConvergenceTestFn(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx);
336: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
337: SLEPC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
338: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedAbsolute;
339: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedRelative;
340: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedNorm;
341: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedMaxIt;
342: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,SVDConvergenceTestFn*,void*,PetscCtxDestroyFn*);
344: /*S
345: SVDStoppingTestFn - A prototype of an SVD stopping test function that would be passed to SVDSetStoppingTestFunction()
347: Calling Sequence:
348: + svd - the singular value solver context
349: . its - current number of iterations
350: . max_it - maximum number of iterations
351: . nconv - number of currently converged singular triplets
352: . nsv - number of requested singular triplets
353: . reason - [output] result of the stopping test
354: - ctx - [optional] user-defined context for private data for the
355: stopping test routine (may be NULL)
357: Level: advanced
359: .seealso: [](ch:svd), `SVDSetStoppingTestFunction()`
360: S*/
361: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDStoppingTestFn(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx);
363: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
364: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
365: SLEPC_EXTERN SVDStoppingTestFn SVDStoppingBasic;
366: SLEPC_EXTERN SVDStoppingTestFn SVDStoppingThreshold;
367: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,SVDStoppingTestFn*,void*,PetscCtxDestroyFn*);
369: /* --------- options specific to particular solvers -------- */
371: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
372: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
373: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
374: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
376: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
377: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
378: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
379: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
381: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
382: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
384: /*E
385: SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
386: TRLanczos GSVD solver
388: Level: advanced
390: .seealso: [](ch:svd), `SVDTRLanczosSetGBidiag()`, `SVDTRLanczosGetGBidiag()`
391: E*/
392: typedef enum {
393: SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
394: SVD_TRLANCZOS_GBIDIAG_UPPER, /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
395: SVD_TRLANCZOS_GBIDIAG_LOWER /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
396: } SVDTRLanczosGBidiag;
397: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];
399: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
400: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
401: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
402: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
403: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
404: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
405: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
406: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
407: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
408: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
409: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
410: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
411: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
412: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);
414: /*E
415: SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library
417: Level: advanced
419: .seealso: [](ch:svd), `SVDPRIMMESetMethod()`, `SVDPRIMMEGetMethod()`
420: E*/
421: typedef enum { SVD_PRIMME_HYBRID = 1,
422: SVD_PRIMME_NORMALEQUATIONS = 2,
423: SVD_PRIMME_AUGMENTED = 3 } SVDPRIMMEMethod;
424: SLEPC_EXTERN const char *SVDPRIMMEMethods[];
426: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
427: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
428: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
429: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);
431: /*E
432: SVDKSVDEigenMethod - determines the method to solve the eigenproblem within KSVD
434: Level: advanced
436: .seealso: [](ch:svd), `SVDKSVDSetEigenMethod()`, `SVDKSVDGetEigenMethod()`
437: E*/
438: typedef enum { SVD_KSVD_EIGEN_MRRR = 1,
439: SVD_KSVD_EIGEN_DC = 2,
440: SVD_KSVD_EIGEN_ELPA = 3 } SVDKSVDEigenMethod;
441: SLEPC_EXTERN const char *SVDKSVDEigenMethods[];
443: /*E
444: SVDKSVDPolarMethod - determines the method to compute the polar decomposition within KSVD
446: Level: advanced
448: .seealso: [](ch:svd), `SVDKSVDSetPolarMethod()`, `SVDKSVDGetPolarMethod()`
449: E*/
450: typedef enum { SVD_KSVD_POLAR_QDWH = 1,
451: SVD_KSVD_POLAR_ZOLOPD = 2 } SVDKSVDPolarMethod;
452: SLEPC_EXTERN const char *SVDKSVDPolarMethods[];
454: SLEPC_EXTERN PetscErrorCode SVDKSVDSetEigenMethod(SVD,SVDKSVDEigenMethod);
455: SLEPC_EXTERN PetscErrorCode SVDKSVDGetEigenMethod(SVD,SVDKSVDEigenMethod*);
456: SLEPC_EXTERN PetscErrorCode SVDKSVDSetPolarMethod(SVD,SVDKSVDPolarMethod);
457: SLEPC_EXTERN PetscErrorCode SVDKSVDGetPolarMethod(SVD,SVDKSVDPolarMethod*);