Actual source code: slepcsvd.h

slepc-main 2024-11-15
Report Typos and Errors
  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,   /* 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_DIVERGED_SYMMETRY_LOST       = -3,
131:               SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;
132: SLEPC_EXTERN const char *const*SVDConvergedReasons;

134: SLEPC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
135: SLEPC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
136: SLEPC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
137: SLEPC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
138: SLEPC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
139: SLEPC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
140: SLEPC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
141: SLEPC_EXTERN PetscErrorCode SVDSetProblemType(SVD,SVDProblemType);
142: SLEPC_EXTERN PetscErrorCode SVDGetProblemType(SVD,SVDProblemType*);
143: SLEPC_EXTERN PetscErrorCode SVDIsGeneralized(SVD,PetscBool*);
144: SLEPC_EXTERN PetscErrorCode SVDIsHyperbolic(SVD,PetscBool*);
145: SLEPC_EXTERN PetscErrorCode SVDSetOperators(SVD,Mat,Mat);
146: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDSetOperators()", ) static inline PetscErrorCode SVDSetOperator(SVD svd,Mat A) {return SVDSetOperators(svd,A,PETSC_NULLPTR);}
147: SLEPC_EXTERN PetscErrorCode SVDGetOperators(SVD,Mat*,Mat*);
148: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDGetOperators()", ) static inline PetscErrorCode SVDGetOperator(SVD svd,Mat *A) {return SVDGetOperators(svd,A,PETSC_NULLPTR);}
149: SLEPC_EXTERN PetscErrorCode SVDSetSignature(SVD,Vec);
150: SLEPC_EXTERN PetscErrorCode SVDGetSignature(SVD,Vec);
151: SLEPC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec[],PetscInt,Vec[]);
152: 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);}
153: 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);}
154: SLEPC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
155: SLEPC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
156: SLEPC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
157: SLEPC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
158: SLEPC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
159: SLEPC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
160: SLEPC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
161: SLEPC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
162: SLEPC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
163: SLEPC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
164: SLEPC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
165: SLEPC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
166: SLEPC_EXTERN PetscErrorCode SVDSetDSType(SVD);
167: SLEPC_EXTERN PetscErrorCode SVDSetUp(SVD);
168: SLEPC_EXTERN PetscErrorCode SVDSolve(SVD);
169: SLEPC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
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 SVDSetStoppingTest(SVD,SVDStop);
177: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
178: SLEPC_EXTERN PetscErrorCode SVDStoppingBasic(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
179: SLEPC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
180: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
181: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
182: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
183: 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);}
184: 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);}
185: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
186: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
187: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
188: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "SVDErrorView()", ) static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
189: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
190: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
191: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
192: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonView()", ) static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
193: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
194: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
195: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
196: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
197: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
198: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
199: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
200: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
201: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
202: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);

204: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
205: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscCtxDestroyFn*);
206: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
207: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);

209: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
210: SLEPC_EXTERN PetscErrorCode SVDMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
211: SLEPC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
212: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
213: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
214: SLEPC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
215: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
216: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
217: SLEPC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
218: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
219: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
220: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
221: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat**);
222: SLEPC_EXTERN PetscErrorCode SVDMonitorConditioning(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);

224: SLEPC_EXTERN PetscFunctionList SVDList;
225: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
226: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
227: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
228: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
229: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));

231: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);

233: /*S
234:   SVDConvergenceTestFn - A prototype of an SVD convergence test function that would be passed to SVDSetConvergenceTestFunction()

236:   Calling Sequence:
237: +   svd    - singular value solver context obtained from SVDCreate()
238: .   sigma  - computed singular value
239: .   res    - residual norm associated to the singular triplet
240: .   errest - [output] computed error estimate
241: -   ctx    - [optional] user-defined context for private data for the
242:              convergence test routine (may be NULL)

244:   Level: advanced

246: .seealso: SVDSetConvergenceTestFunction()
247: S*/
248: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SVDConvergenceTestFn)(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx);

250: /*S
251:   SVDStoppingTestFn - A prototype of an SVD stopping test function that would be passed to SVDSetStoppingTestFunction()

253:   Calling Sequence:
254: +   svd    - singular value solver context obtained from SVDCreate()
255: .   its    - current number of iterations
256: .   max_it - maximum number of iterations
257: .   nconv  - number of currently converged singular triplets
258: .   nsv    - number of requested singular triplets
259: .   reason - [output] result of the stopping test
260: -   ctx    - [optional] user-defined context for private data for the
261:              stopping test routine (may be NULL)

263:   Level: advanced

265: .seealso: SVDSetStoppingTestFunction()
266: S*/
267: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SVDStoppingTestFn)(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx);

269: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,SVDConvergenceTestFn*,void*,PetscCtxDestroyFn*);
270: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,SVDStoppingTestFn*,void*,PetscCtxDestroyFn*);

272: /* --------- options specific to particular solvers -------- */

274: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
275: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
276: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
277: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);

279: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
280: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
281: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
282: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);

284: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
285: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

287: /*E
288:     SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
289:     TRLanczos GSVD solver

291:     Level: advanced

293: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
294: E*/
295: typedef enum {
296:   SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
297:   SVD_TRLANCZOS_GBIDIAG_UPPER,  /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
298:   SVD_TRLANCZOS_GBIDIAG_LOWER   /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
299: } SVDTRLanczosGBidiag;
300: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];

302: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
303: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
304: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
305: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
306: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
307: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
308: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
309: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
310: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
311: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
312: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
313: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
314: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
315: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);

317: /*E
318:     SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library

320:     Level: advanced

322: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
323: E*/
324: typedef enum { SVD_PRIMME_HYBRID=1,
325:                SVD_PRIMME_NORMALEQUATIONS,
326:                SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
327: SLEPC_EXTERN const char *SVDPRIMMEMethods[];

329: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
330: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
331: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
332: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);

334: /*E
335:     SVDKSVDEigenMethod - determines the method to solve the eigenproblem within KSVD

337:     Level: advanced

339: .seealso: SVDKSVDSetEigenMethod(), SVDKSVDGetEigenMethod()
340: E*/
341: typedef enum { SVD_KSVD_EIGEN_MRRR=1,
342:                SVD_KSVD_EIGEN_DC,
343:                SVD_KSVD_EIGEN_ELPA } SVDKSVDEigenMethod;
344: SLEPC_EXTERN const char *SVDKSVDEigenMethods[];

346: /*E
347:     SVDKSVDPolarMethod - determines the method to compute the polar decomposition within KSVD

349:     Level: advanced

351: .seealso: SVDKSVDSetPolarMethod(), SVDKSVDGetPolarMethod()
352: E*/
353: typedef enum { SVD_KSVD_POLAR_QDWH=1,
354:                SVD_KSVD_POLAR_ZOLOPD } SVDKSVDPolarMethod;
355: SLEPC_EXTERN const char *SVDKSVDPolarMethods[];

357: SLEPC_EXTERN PetscErrorCode SVDKSVDSetEigenMethod(SVD,SVDKSVDEigenMethod);
358: SLEPC_EXTERN PetscErrorCode SVDKSVDGetEigenMethod(SVD,SVDKSVDEigenMethod*);
359: SLEPC_EXTERN PetscErrorCode SVDKSVDSetPolarMethod(SVD,SVDKSVDPolarMethod);
360: SLEPC_EXTERN PetscErrorCode SVDKSVDGetPolarMethod(SVD,SVDKSVDPolarMethod*);