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 - 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 SVDGetConvergedReason(SVD,SVDConvergedReason*);
191: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
192: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
193: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
194: 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);}
195: 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);}
196: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
197: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
198: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
199: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "SVDErrorView()", ) static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
200: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
201: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
202: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
203: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonView()", ) static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
204: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
205: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
206: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
207: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
208: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
209: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
210: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
211: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
212: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
213: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
215: /*S
216: SVDMonitorFn - A function prototype for functions provided to SVDMonitorSet()
218: Calling Sequence:
219: + svd - singular value solver context obtained from SVDCreate()
220: . its - iteration number
221: . nconv - number of converged singular triplets
222: . sigma - singular values
223: . errest - relative error estimates for each singular triplet
224: . nest - number of error estimates
225: - ctx - optional monitoring context, as provided with SVDMonitorSet()
227: Level: beginner
229: .seealso: SVDMonitorSet()
230: S*/
231: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorFn(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *ctx);
233: /*S
234: SVDMonitorRegisterFn - A function prototype for functions provided to SVDMonitorRegister()
236: Calling Sequence:
237: + svd - singular value solver context obtained from SVDCreate()
238: . its - iteration number
239: . nconv - number of converged singular triplets
240: . sigma - singular values
241: . errest - relative error estimates for each singular triplet
242: . nest - number of error estimates
243: - ctx - PetscViewerAndFormat object
245: Level: beginner
247: Note:
248: This is an SVDMonitorFn specialized for a context of PetscViewerAndFormat.
250: .seealso: SVDMonitorSet(), SVDMonitorRegister(), SVDMonitorFn, SVDMonitorRegisterCreateFn, SVDMonitorRegisterDestroyFn
251: S*/
252: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterFn(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *ctx);
254: /*S
255: SVDMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to SVDMonitorRegister()
257: Calling Sequence:
258: + viewer - the viewer to be used with the SVDMonitorRegisterFn
259: . format - the format of the viewer
260: . ctx - a context for the monitor
261: - result - a PetscViewerAndFormat object
263: Level: beginner
265: .seealso: SVDMonitorRegisterFn, SVDMonitorSet(), SVDMonitorRegister(), SVDMonitorFn, SVDMonitorRegisterDestroyFn
266: S*/
267: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **result);
269: /*S
270: SVDMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to SVDMonitorRegister()
272: Calling Sequence:
273: . vf - a PetscViewerAndFormat object to be destroyed, including any context
275: Level: beginner
277: .seealso: SVDMonitorRegisterFn, SVDMonitorSet(), SVDMonitorRegister(), SVDMonitorFn, SVDMonitorRegisterCreateFn
278: S*/
279: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterDestroyFn(PetscViewerAndFormat **result);
281: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
282: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,SVDMonitorFn,void*,PetscCtxDestroyFn*);
283: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
284: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);
286: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
287: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorFirst;
288: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorFirstDrawLG;
289: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorFirstDrawLGCreate;
290: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorAll;
291: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorAllDrawLG;
292: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorAllDrawLGCreate;
293: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorConverged;
294: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorConvergedCreate;
295: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorConvergedDrawLG;
296: SLEPC_EXTERN SVDMonitorRegisterCreateFn SVDMonitorConvergedDrawLGCreate;
297: SLEPC_EXTERN SVDMonitorRegisterDestroyFn SVDMonitorConvergedDestroy;
298: SLEPC_EXTERN SVDMonitorRegisterFn SVDMonitorConditioning;
300: SLEPC_EXTERN PetscFunctionList SVDList;
301: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
302: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
303: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
304: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
305: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,SVDMonitorRegisterFn*,SVDMonitorRegisterCreateFn*,SVDMonitorRegisterDestroyFn*);
307: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
308: SLEPC_EXTERN PetscErrorCode SVDReallocateSolution(SVD,PetscInt);
310: /*S
311: SVDConvergenceTestFn - A prototype of an SVD convergence test function that would be passed to SVDSetConvergenceTestFunction()
313: Calling Sequence:
314: + svd - singular value solver context obtained from SVDCreate()
315: . sigma - computed singular value
316: . res - residual norm associated to the singular triplet
317: . errest - [output] computed error estimate
318: - ctx - [optional] user-defined context for private data for the
319: convergence test routine (may be NULL)
321: Level: advanced
323: .seealso: SVDSetConvergenceTestFunction()
324: S*/
325: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDConvergenceTestFn(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx);
327: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
328: SLEPC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
329: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedAbsolute;
330: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedRelative;
331: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedNorm;
332: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedMaxIt;
333: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,SVDConvergenceTestFn*,void*,PetscCtxDestroyFn*);
335: /*S
336: SVDStoppingTestFn - A prototype of an SVD stopping test function that would be passed to SVDSetStoppingTestFunction()
338: Calling Sequence:
339: + svd - singular value solver context obtained from SVDCreate()
340: . its - current number of iterations
341: . max_it - maximum number of iterations
342: . nconv - number of currently converged singular triplets
343: . nsv - number of requested singular triplets
344: . reason - [output] result of the stopping test
345: - ctx - [optional] user-defined context for private data for the
346: stopping test routine (may be NULL)
348: Level: advanced
350: .seealso: SVDSetStoppingTestFunction()
351: S*/
352: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDStoppingTestFn(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx);
354: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
355: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
356: SLEPC_EXTERN SVDStoppingTestFn SVDStoppingBasic;
357: SLEPC_EXTERN SVDStoppingTestFn SVDStoppingThreshold;
358: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,SVDStoppingTestFn*,void*,PetscCtxDestroyFn*);
360: /* --------- options specific to particular solvers -------- */
362: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
363: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
364: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
365: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
367: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
368: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
369: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
370: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
372: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
373: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
375: /*E
376: SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
377: TRLanczos GSVD solver
379: Level: advanced
381: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
382: E*/
383: typedef enum {
384: SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
385: SVD_TRLANCZOS_GBIDIAG_UPPER, /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
386: SVD_TRLANCZOS_GBIDIAG_LOWER /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
387: } SVDTRLanczosGBidiag;
388: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];
390: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
391: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
392: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
393: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
394: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
395: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
396: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
397: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
398: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
399: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
400: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
401: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
402: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
403: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);
405: /*E
406: SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library
408: Level: advanced
410: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
411: E*/
412: typedef enum { SVD_PRIMME_HYBRID = 1,
413: SVD_PRIMME_NORMALEQUATIONS = 2,
414: SVD_PRIMME_AUGMENTED = 3 } SVDPRIMMEMethod;
415: SLEPC_EXTERN const char *SVDPRIMMEMethods[];
417: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
418: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
419: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
420: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);
422: /*E
423: SVDKSVDEigenMethod - determines the method to solve the eigenproblem within KSVD
425: Level: advanced
427: .seealso: SVDKSVDSetEigenMethod(), SVDKSVDGetEigenMethod()
428: E*/
429: typedef enum { SVD_KSVD_EIGEN_MRRR = 1,
430: SVD_KSVD_EIGEN_DC = 2,
431: SVD_KSVD_EIGEN_ELPA = 3 } SVDKSVDEigenMethod;
432: SLEPC_EXTERN const char *SVDKSVDEigenMethods[];
434: /*E
435: SVDKSVDPolarMethod - determines the method to compute the polar decomposition within KSVD
437: Level: advanced
439: .seealso: SVDKSVDSetPolarMethod(), SVDKSVDGetPolarMethod()
440: E*/
441: typedef enum { SVD_KSVD_POLAR_QDWH = 1,
442: SVD_KSVD_POLAR_ZOLOPD = 2 } SVDKSVDPolarMethod;
443: SLEPC_EXTERN const char *SVDKSVDPolarMethods[];
445: SLEPC_EXTERN PetscErrorCode SVDKSVDSetEigenMethod(SVD,SVDKSVDEigenMethod);
446: SLEPC_EXTERN PetscErrorCode SVDKSVDGetEigenMethod(SVD,SVDKSVDEigenMethod*);
447: SLEPC_EXTERN PetscErrorCode SVDKSVDSetPolarMethod(SVD,SVDKSVDPolarMethod);
448: SLEPC_EXTERN PetscErrorCode SVDKSVDGetPolarMethod(SVD,SVDKSVDPolarMethod*);