Actual source code: slepcsvd.h

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 } SVDErrorType;
 84: SLEPC_EXTERN const char *SVDErrorTypes[];

 86: /*E
 87:     SVDConv - Determines the convergence test

 89:     Level: intermediate

 91: .seealso: SVDSetConvergenceTest(), SVDSetConvergenceTestFunction()
 92: E*/
 93: typedef enum { SVD_CONV_ABS,
 94:                SVD_CONV_REL,
 95:                SVD_CONV_NORM,
 96:                SVD_CONV_MAXIT,
 97:                SVD_CONV_USER } SVDConv;

 99: /*E
100:     SVDStop - Determines the stopping test

102:     Level: advanced

104: .seealso: SVDSetStoppingTest(), SVDSetStoppingTestFunction()
105: E*/
106: typedef enum { SVD_STOP_BASIC,
107:                SVD_STOP_USER } SVDStop;

109: /*E
110:     SVDConvergedReason - Reason a singular value solver was said to
111:          have converged or diverged

113:    Level: intermediate

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

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

195: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
196: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
197: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
198: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);

200: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
201: SLEPC_EXTERN PetscErrorCode SVDMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
202: SLEPC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
203: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
204: SLEPC_EXTERN PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
205: SLEPC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
206: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
207: SLEPC_EXTERN PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
208: SLEPC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
209: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
210: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
211: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
212: SLEPC_EXTERN PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat**);

214: SLEPC_EXTERN PetscFunctionList SVDList;
215: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
216: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
217: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
218: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
219: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));

221: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);

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

225: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
226: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
227: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
228: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);

230: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
231: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
232: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
233: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);

235: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
236: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

238: /*E
239:     SVDTRLanczosGBidiag - determines the bidiagonalization choice for the
240:     TRLanczos GSVD solver

242:     Level: advanced

244: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGetGBidiag()
245: E*/
246: typedef enum {
247:   SVD_TRLANCZOS_GBIDIAG_SINGLE, /* single bidiagonalization (Qa) */
248:   SVD_TRLANCZOS_GBIDIAG_UPPER,  /* joint bidiagonalization, both Qa and Qb in upper bidiagonal form */
249:   SVD_TRLANCZOS_GBIDIAG_LOWER   /* joint bidiagonalization, Qa lower bidiagonal, Qb upper bidiagonal */
250: } SVDTRLanczosGBidiag;
251: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];

253: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
254: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
255: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
256: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
257: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
258: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
259: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
260: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
261: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
262: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
263: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
264: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);

266: /*E
267:     SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library

269:     Level: advanced

271: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
272: E*/
273: typedef enum { SVD_PRIMME_HYBRID=1,
274:                SVD_PRIMME_NORMALEQUATIONS,
275:                SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
276: SLEPC_EXTERN const char *SVDPRIMMEMethods[];

278: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
279: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
280: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
281: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);

283: #endif