Actual source code: slepcsvd.h
slepc-3.17.1 2022-04-11
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
16: #include <slepceps.h>
17: #include <slepcbv.h>
18: #include <slepcds.h>
20: SLEPC_EXTERN PetscErrorCode SVDInitializePackage(void);
22: /*S
23: SVD - Abstract SLEPc object that manages all the singular value
24: problem solvers.
26: Level: beginner
28: .seealso: SVDCreate()
29: S*/
30: typedef struct _p_SVD* SVD;
32: /*J
33: SVDType - String with the name of a SLEPc singular value solver
35: Level: beginner
37: .seealso: SVDSetType(), SVD
38: J*/
39: typedef const char* SVDType;
40: #define SVDCROSS "cross"
41: #define SVDCYCLIC "cyclic"
42: #define SVDLAPACK "lapack"
43: #define SVDLANCZOS "lanczos"
44: #define SVDTRLANCZOS "trlanczos"
45: #define SVDRANDOMIZED "randomized"
46: #define SVDSCALAPACK "scalapack"
47: #define SVDELEMENTAL "elemental"
48: #define SVDPRIMME "primme"
50: /* Logging support */
51: SLEPC_EXTERN PetscClassId SVD_CLASSID;
53: /*E
54: SVDProblemType - Determines the type of singular value problem
56: Level: beginner
58: .seealso: SVDSetProblemType(), SVDGetProblemType()
59: E*/
60: typedef enum { SVD_STANDARD=1,
61: SVD_GENERALIZED /* GSVD */
62: } SVDProblemType;
64: /*E
65: SVDWhich - Determines whether largest or smallest singular triplets
66: are to be computed
68: Level: intermediate
70: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
71: E*/
72: typedef enum { SVD_LARGEST,
73: SVD_SMALLEST } SVDWhich;
75: /*E
76: SVDErrorType - The error type used to assess accuracy of computed solutions
78: Level: intermediate
80: .seealso: SVDComputeError()
81: E*/
82: typedef enum { SVD_ERROR_ABSOLUTE,
83: SVD_ERROR_RELATIVE,
84: SVD_ERROR_NORM } SVDErrorType;
85: SLEPC_EXTERN const char *SVDErrorTypes[];
87: /*E
88: SVDConv - Determines the convergence test
90: Level: intermediate
92: .seealso: SVDSetConvergenceTest(), SVDSetConvergenceTestFunction()
93: E*/
94: typedef enum { SVD_CONV_ABS,
95: SVD_CONV_REL,
96: SVD_CONV_NORM,
97: SVD_CONV_MAXIT,
98: SVD_CONV_USER } SVDConv;
100: /*E
101: SVDStop - Determines the stopping test
103: Level: advanced
105: .seealso: SVDSetStoppingTest(), SVDSetStoppingTestFunction()
106: E*/
107: typedef enum { SVD_STOP_BASIC,
108: SVD_STOP_USER } SVDStop;
110: /*E
111: SVDConvergedReason - Reason a singular value solver was said to
112: have converged or diverged
114: Level: intermediate
116: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
117: E*/
118: typedef enum {/* converged */
119: SVD_CONVERGED_TOL = 1,
120: SVD_CONVERGED_USER = 2,
121: SVD_CONVERGED_MAXIT = 3,
122: /* diverged */
123: SVD_DIVERGED_ITS = -1,
124: SVD_DIVERGED_BREAKDOWN = -2,
125: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
126: SLEPC_EXTERN const char *const*SVDConvergedReasons;
128: SLEPC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
129: SLEPC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
130: SLEPC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
131: SLEPC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
132: SLEPC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
133: SLEPC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
134: SLEPC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
135: SLEPC_EXTERN PetscErrorCode SVDSetProblemType(SVD,SVDProblemType);
136: SLEPC_EXTERN PetscErrorCode SVDGetProblemType(SVD,SVDProblemType*);
137: SLEPC_EXTERN PetscErrorCode SVDIsGeneralized(SVD,PetscBool*);
138: SLEPC_EXTERN PetscErrorCode SVDSetOperators(SVD,Mat,Mat);
139: PETSC_DEPRECATED_FUNCTION("Use SVDSetOperators()") static inline PetscErrorCode SVDSetOperator(SVD svd,Mat A) {return SVDSetOperators(svd,A,NULL);}
140: SLEPC_EXTERN PetscErrorCode SVDGetOperators(SVD,Mat*,Mat*);
141: PETSC_DEPRECATED_FUNCTION("Use SVDGetOperators()") static inline PetscErrorCode SVDGetOperator(SVD svd,Mat *A) {return SVDGetOperators(svd,A,NULL);}
142: SLEPC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec[],PetscInt,Vec[]);
143: PETSC_DEPRECATED_FUNCTION("Use SVDSetInitialSpaces()") static inline PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt nr,Vec *isr) {return SVDSetInitialSpaces(svd,nr,isr,0,NULL);}
144: PETSC_DEPRECATED_FUNCTION("Use SVDSetInitialSpaces()") static inline PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt nl,Vec *isl) {return SVDSetInitialSpaces(svd,0,NULL,nl,isl);}
145: SLEPC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
146: SLEPC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
147: SLEPC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
148: SLEPC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
149: SLEPC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
150: SLEPC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
151: SLEPC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
152: SLEPC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
153: SLEPC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
154: SLEPC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
155: SLEPC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
156: SLEPC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
157: SLEPC_EXTERN PetscErrorCode SVDSetUp(SVD);
158: SLEPC_EXTERN PetscErrorCode SVDSolve(SVD);
159: SLEPC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
160: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,PetscErrorCode (*)(SVD,PetscReal,PetscReal,PetscReal*,void*),void*,PetscErrorCode (*)(void*));
161: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
162: SLEPC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
163: SLEPC_EXTERN PetscErrorCode SVDConvergedAbsolute(SVD,PetscReal,PetscReal,PetscReal*,void*);
164: SLEPC_EXTERN PetscErrorCode SVDConvergedRelative(SVD,PetscReal,PetscReal,PetscReal*,void*);
165: SLEPC_EXTERN PetscErrorCode SVDConvergedNorm(SVD,PetscReal,PetscReal,PetscReal*,void*);
166: SLEPC_EXTERN PetscErrorCode SVDConvergedMaxIt(SVD,PetscReal,PetscReal,PetscReal*,void*);
167: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
168: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
169: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
170: SLEPC_EXTERN PetscErrorCode SVDStoppingBasic(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
171: SLEPC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
172: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
173: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
174: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
175: PETSC_DEPRECATED_FUNCTION("Use SVDComputeError()") static inline PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *r) {return SVDComputeError(svd,i,SVD_ERROR_RELATIVE,r);}
176: 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);}
177: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
178: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
179: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
180: PETSC_DEPRECATED_FUNCTION("Use SVDErrorView()") static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
181: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
182: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
183: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
184: PETSC_DEPRECATED_FUNCTION("Use SVDConvergedReasonView() (since version 3.14)") static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
185: PETSC_DEPRECATED_FUNCTION("Use SVDConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
186: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
187: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
188: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
189: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
190: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
191: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
192: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
193: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
194: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
196: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
197: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
198: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
199: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);
201: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
202: SLEPC_EXTERN PetscErrorCode SVDMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
203: SLEPC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
204: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
205: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
206: SLEPC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
207: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
208: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
209: SLEPC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
210: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
211: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
212: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
213: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat**);
215: SLEPC_EXTERN PetscFunctionList SVDList;
216: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
217: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
218: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
219: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
220: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));
222: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
224: /* --------- options specific to particular solvers -------- */
226: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
227: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
228: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
229: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
231: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
232: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
233: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
234: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
236: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
237: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
239: /*E
240: SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
241: TRLanczos GSVD solver
243: Level: advanced
245: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
246: E*/
247: typedef enum {
248: SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
249: SVD_TRLANCZOS_GBIDIAG_UPPER, /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
250: SVD_TRLANCZOS_GBIDIAG_LOWER /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
251: } SVDTRLanczosGBidiag;
252: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];
254: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
255: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
256: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
257: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
258: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
259: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
260: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
261: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
262: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
263: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
264: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
265: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
267: /*E
268: SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library
270: Level: advanced
272: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
273: E*/
274: typedef enum { SVD_PRIMME_HYBRID=1,
275: SVD_PRIMME_NORMALEQUATIONS,
276: SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
277: SLEPC_EXTERN const char *SVDPRIMMEMethods[];
279: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
280: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
281: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
282: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);
284: #endif