Actual source code: slepcsvd.h

slepc-3.18.2 2023-01-26
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: #if !defined(SLEPCSVD_H)
 15: #define SLEPCSVD_H

 17: #include <slepceps.h>
 18: #include <slepcbv.h>
 19: #include <slepcds.h>

 21: /* SUBMANSEC = SVD */

 23: SLEPC_EXTERN PetscErrorCode SVDInitializePackage(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 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: SVDSetProblemType(), SVDGetProblemType()
 62: E*/
 63: typedef enum { SVD_STANDARD=1,
 64:                SVD_GENERALIZED,   /* GSVD */
 65:                SVD_HYPERBOLIC     /* HSVD */
 66:              } SVDProblemType;

 68: /*E
 69:     SVDWhich - Determines whether largest or smallest singular triplets
 70:     are to be computed

 72:     Level: intermediate

 74: .seealso: 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: 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: 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: SVDSetStoppingTest(), SVDSetStoppingTestFunction()
110: E*/
111: typedef enum { SVD_STOP_BASIC,
112:                SVD_STOP_USER } SVDStop;

114: /*E
115:     SVDConvergedReason - Reason a singular value solver was said to
116:          have converged or diverged

118:    Level: intermediate

120: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
121: E*/
122: typedef enum {/* converged */
123:               SVD_CONVERGED_TOL                =  1,
124:               SVD_CONVERGED_USER               =  2,
125:               SVD_CONVERGED_MAXIT              =  3,
126:               /* diverged */
127:               SVD_DIVERGED_ITS                 = -1,
128:               SVD_DIVERGED_BREAKDOWN           = -2,
129:               SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;
130: SLEPC_EXTERN const char *const*SVDConvergedReasons;

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

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

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

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

230: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);

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

234: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
235: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
236: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
237: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);

239: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
240: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
241: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
242: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);

244: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
245: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

247: /*E
248:     SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
249:     TRLanczos GSVD solver

251:     Level: advanced

253: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
254: E*/
255: typedef enum {
256:   SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
257:   SVD_TRLANCZOS_GBIDIAG_UPPER,  /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
258:   SVD_TRLANCZOS_GBIDIAG_LOWER   /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
259: } SVDTRLanczosGBidiag;
260: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];

262: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
263: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
264: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
265: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
266: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
267: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
268: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
269: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
270: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
271: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
272: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
273: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
274: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
275: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);

277: /*E
278:     SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library

280:     Level: advanced

282: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
283: E*/
284: typedef enum { SVD_PRIMME_HYBRID=1,
285:                SVD_PRIMME_NORMALEQUATIONS,
286:                SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
287: SLEPC_EXTERN const char *SVDPRIMMEMethods[];

289: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
290: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
291: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
292: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);

294: #endif