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 - SLEPc object that manages all the singular value problem solvers.

 28:    Level: beginner

 30: .seealso: [](ch:svd), `SVDCreate()`
 31: S*/
 32: typedef struct _p_SVD* SVD;

 34: /*J
 35:    SVDType - String with the name of a singular value solver.

 37:    Level: beginner

 39: .seealso: [](ch:svd), `SVDSetType()`, `SVD`
 40: J*/
 41: typedef const char *SVDType;
 42: #define SVDCROSS       "cross"
 43: #define SVDCYCLIC      "cyclic"
 44: #define SVDLANCZOS     "lanczos"
 45: #define SVDTRLANCZOS   "trlanczos"
 46: #define SVDRANDOMIZED  "randomized"
 47: #define SVDLAPACK      "lapack"
 48: #define SVDSCALAPACK   "scalapack"
 49: #define SVDKSVD        "ksvd"
 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 the singular value problem.

 59:    Values:
 60: +  `SVD_STANDARD`    - standard SVD (SVD)
 61: .  `SVD_GENERALIZED` - generalized SVD (GSVD)
 62: -  `SVD_HYPERBOLIC`  - hyperbolic SVD (HSVD)

 64:    Level: beginner

 66: .seealso: [](ch:svd), `SVDSetProblemType()`, `SVDGetProblemType()`
 67: E*/
 68: typedef enum { SVD_STANDARD    = 1,
 69:                SVD_GENERALIZED = 2,
 70:                SVD_HYPERBOLIC  = 3
 71:              } SVDProblemType;

 73: /*MC
 74:    SVD_STANDARD - A (standard) singular value problem (SVD).

 76:    Note:
 77:    The problem is formulated as $A=U\Sigma V^*$, where $A$ is a possibly
 78:    non-square matrix.

 80:    Level: beginner

 82: .seealso: [](ch:svd), `SVDProblemType`, `SVDSetProblemType()`, `SVD_GENERALIZED`, `SVD_HYPERBOLIC`
 83: M*/

 85: /*MC
 86:    SVD_GENERALIZED - A generalized singular value problem (GSVD).

 88:    Note:
 89:    The problem is formulated as $U^*AX=C$, $V^*BX=S$, where $A$ and $B$
 90:    have the same number of columns.

 92:    Level: beginner

 94: .seealso: [](ch:svd), `SVDProblemType`, `SVDSetProblemType()`, `SVD_STANDARD`, `SVD_HYPERBOLIC`
 95: M*/

 97: /*MC
 98:    SVD_HYPERBOLIC - A hyperbolic singular value problem (HSVD).

100:    Note:
101:    The problem is formulated as $A=U\Sigma V^*$, with $U^*\Omega U=\tilde\Omega$,
102:    where $A$ is a possibly non-square matrix, and $\Omega$, $\tilde\Omega$
103:    are signature matrices.

105:    Level: beginner

107: .seealso: [](ch:svd), `SVDProblemType`, `SVDSetProblemType()`, `SVD_STANDARD`, `SVD_GENERALIZED`
108: M*/

110: /*E
111:    SVDWhich - Determines whether largest or smallest singular values
112:    are to be computed.

114:    Values:
115: +  `SVD_LARGEST`  - largest singular values
116: -  `SVD_SMALLEST` - smallest singular values

118:    Level: intermediate

120: .seealso: [](ch:svd), `SVDSetWhichSingularTriplets()`, `SVDGetWhichSingularTriplets()`
121: E*/
122: typedef enum { SVD_LARGEST,
123:                SVD_SMALLEST } SVDWhich;

125: /*MC
126:    SVD_LARGEST - The solver is configured to compute largest singular values.

128:    Note:
129:    This is the default.

131:    Level: intermediate

133: .seealso: [](ch:svd), `SVDWhich`, `SVDSetWhichSingularTriplets()`, `SVD_SMALLEST`
134: M*/

136: /*MC
137:    SVD_SMALLEST - The solver is configured to compute smallest singular values.

139:    Note:
140:    Computing small singular values is generally more difficult than computing
141:    largest ones, because in many cases these values are very small and
142:    tightly clustered together. In the case of rank-deficient matrices, smallest
143:    singular values are zero, and this may pose difficulties to the solvers.

145:    Level: intermediate

147: .seealso: [](ch:svd), `SVDWhich`, `SVDSetWhichSingularTriplets()`, `SVD_LARGEST`
148: M*/

150: /*E
151:    SVDErrorType - The error type used to assess accuracy of computed solutions.

153:    Values:
154: +  `SVD_ERROR_ABSOLUTE` - compute error bound as $\|r\|$
155: .  `SVD_ERROR_RELATIVE` - compute error bound as $\|r\|/\sigma$
156: -  `SVD_ERROR_NORM`     - compute error bound as $\|r\|/\max\{\|A\|,\|B\|\}$

158:    Note:
159:    The residual norm $\|r\|$ is actually computed from two parts, such as
160:    $\sqrt{\eta_1^2+\eta_2^2}$ with $\eta_1 = \|Av-\sigma u\|_2$ and
161:    $\eta_2 = \|A^*u-\sigma v\|_2$, see more details at `SVDComputeError()`.
162:    There is also a normalization factor related to the norm of the vectors,
163:    which also varies with the problem type.

165:    Level: intermediate

167: .seealso: [](ch:svd), `SVDComputeError()`, `SVDProblemType`
168: E*/
169: typedef enum { SVD_ERROR_ABSOLUTE,
170:                SVD_ERROR_RELATIVE,
171:                SVD_ERROR_NORM } SVDErrorType;
172: SLEPC_EXTERN const char *SVDErrorTypes[];

174: /*E
175:    SVDConv - The convergence criterion to be used by the solver.

177:    Values:
178: +  `SVD_CONV_ABS`   - absolute convergence criterion, $\|r\|$
179: .  `SVD_CONV_REL`   - convergence criterion relative to singular value, $\|r\|/\sigma$
180: .  `SVD_CONV_NORM`  - convergence criterion relative to matrix norms, $\|r\|/\max\{\|A\|,\|B\|\}$
181: .  `SVD_CONV_MAXIT` - no convergence until maximum number of iterations has been reached
182: -  `SVD_CONV_USER`  - convergence dictated by user-provided function

184:    Note:
185:    The `SVD_CONV_MAXIT` convergence criterion is used only in `SVDRANDOMIZED`.

187:    Level: intermediate

189: .seealso: [](ch:svd), `SVDSetConvergenceTest()`, `SVDSetConvergenceTestFunction()`, `SVDSetTolerances()`
190: E*/
191: typedef enum { SVD_CONV_ABS,
192:                SVD_CONV_REL,
193:                SVD_CONV_NORM,
194:                SVD_CONV_MAXIT,
195:                SVD_CONV_USER } SVDConv;

197: /*E
198:    SVDStop - The stopping test to decide the termination of the outer loop
199:    of the singular value solver.

201:    Values:
202: +  `SVD_STOP_BASIC`     - default stopping test
203: .  `SVD_STOP_USER`      - user-provided stopping test
204: -  `SVD_STOP_THRESHOLD` - threshold stopping test

206:    Level: advanced

208: .seealso: [](ch:svd), `SVDSetStoppingTest()`, `SVDSetStoppingTestFunction()`
209: E*/
210: typedef enum { SVD_STOP_BASIC,
211:                SVD_STOP_USER,
212:                SVD_STOP_THRESHOLD } SVDStop;

214: /*MC
215:    SVD_STOP_BASIC - The default stopping test.

217:    Note:
218:    By default, the termination of the outer loop is decided by calling
219:    `SVDStoppingBasic()`, which will stop if all requested singular values are converged,
220:    or if the maximum number of iterations has been reached.

222:    Level: advanced

224: .seealso: [](ch:svd), `SVDStop`, `SVDSetStoppingTest()`, `SVDStoppingBasic()`
225: M*/

227: /*MC
228:    SVD_STOP_USER - The user-provided stopping test.

230:    Note:
231:    Customized stopping test using the user-provided function given with
232:    `SVDSetStoppingTestFunction()`.

234:    Level: advanced

236: .seealso: [](ch:svd), `SVDStop`, `SVDSetStoppingTest()`, `SVDSetStoppingTestFunction()`
237: M*/

239: /*MC
240:    SVD_STOP_THRESHOLD - The threshold stopping test.

242:    Note:
243:    When a threshold has been provided with `SVDSetThreshold()`, the termination
244:    of the outer loop is decided by calling `SVDStoppingThreshold()`, which will
245:    stop when one of the computed singular values is not above/below the threshold.
246:    If a number of wanted singular values has been specified via `SVDSetDimensions()`
247:    then it is also taken into account, and the solver will stop when one of the
248:    two conditions (threshold or number of converged values) is met.

250:    Level: advanced

252: .seealso: [](ch:svd), `SVDStop`, `SVDSetStoppingTest()`, `SVDStoppingThreshold()`, `SVDSetThreshold()`, `SVDSetDimensions()`
253: M*/

255: /*E
256:    SVDConvergedReason - Reason a singular value solver was determined to have
257:    converged or diverged.

259:    Values:
260: +  `SVD_CONVERGED_TOL`          - converged up to tolerance
261: .  `SVD_CONVERGED_USER`         - converged due to a user-defined condition
262: .  `SVD_CONVERGED_MAXIT`        - reached maximum number of iterations with `SVD_CONV_MAXIT` criterion
263: .  `SVD_DIVERGED_ITS`           - exceeded the maximum number of allowed iterations
264: .  `SVD_DIVERGED_BREAKDOWN`     - generic breakdown in method
265: .  `SVD_DIVERGED_SYMMETRY_LOST` - underlying indefinite eigensolver was not able to keep symmetry
266: -  `SVD_CONVERGED_ITERATING`    - the solver is still running

268:    Level: intermediate

270: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDSetTolerances()`
271: E*/
272: typedef enum {/* converged */
273:               SVD_CONVERGED_TOL                =  1,
274:               SVD_CONVERGED_USER               =  2,
275:               SVD_CONVERGED_MAXIT              =  3,
276:               /* diverged */
277:               SVD_DIVERGED_ITS                 = -1,
278:               SVD_DIVERGED_BREAKDOWN           = -2,
279:               SVD_DIVERGED_SYMMETRY_LOST       = -3,
280:               SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;
281: SLEPC_EXTERN const char *const*SVDConvergedReasons;

283: /*MC
284:    SVD_CONVERGED_TOL - The computed error estimates, based on residual norms,
285:    for all requested singular values are below the tolerance.

287:    Level: intermediate

289: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`
290: M*/

292: /*MC
293:    SVD_CONVERGED_USER - The solver was declared converged due to a user-defined condition.

295:    Note:
296:    This happens only when a user-defined stopping test has been set with
297:    `SVDSetStoppingTestFunction()`.

299:    Level: intermediate

301: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`, `SVDSetStoppingTestFunction()`
302: M*/

304: /*MC
305:    SVD_CONVERGED_MAXIT - The solver has reached the maximum number of iterations
306:    with the `SVD_CONV_MAXIT` criterion.

308:    Note:
309:    This is considered a successful exit, because the user wanted to do a fixed
310:    number of iterations. But be aware that the computed solution may be inaccurate,
311:    in particular, individual singular vectors will not have good residual. This
312:    is available in `SVDRANDOMIZED` only.

314:    Level: intermediate

316: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`, `SVD_CONV_MAXIT`, `SVDRANDOMIZED`
317: M*/

319: /*MC
320:    SVD_DIVERGED_ITS - Exceeded the maximum number of allowed iterations
321:    before the convergence criterion was satisfied.

323:    Level: intermediate

325: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`
326: M*/

328: /*MC
329:    SVD_DIVERGED_BREAKDOWN - A breakdown in the solver was detected so the
330:    method could not continue.

332:    Level: intermediate

334: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`
335: M*/

337: /*MC
338:    SVD_DIVERGED_SYMMETRY_LOST - The selected solver uses a pseudo-Lanczos recurrence,
339:    which is numerically unstable, and a symmetry test revealed that instability
340:    had appeared so the solver could not continue.

342:    Level: intermediate

344: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`
345: M*/

347: /*MC
348:    SVD_CONVERGED_ITERATING - This value is returned if `SVDGetConvergedReason()` is called
349:    while `SVDSolve()` is still running.

351:    Level: intermediate

353: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConvergedReason()`, `SVDConvergedReason`
354: M*/

356: /*S
357:    SVDStoppingCtx - Data structure (C struct) to hold additional information to
358:    be used in some stopping test functions.

360:    Level: advanced

362: .seealso: [](ch:svd), `SVDSetStoppingTestFunction()`
363: S*/
364: struct _n_SVDStoppingCtx {
365:   PetscReal firstsv;    /* the value of the first converged singular value */
366:   PetscReal lastsv;     /* the value of the last converged singular value */
367:   PetscReal thres;      /* threshold set with SVDSetThreshold() */
368:   PetscBool threlative; /* threshold is relative */
369:   SVDWhich  which;      /* which singular values are being computed */
370: };
371: typedef struct _n_SVDStoppingCtx* SVDStoppingCtx;

373: SLEPC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
374: SLEPC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
375: SLEPC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
376: SLEPC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
377: SLEPC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
378: SLEPC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
379: SLEPC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
380: SLEPC_EXTERN PetscErrorCode SVDSetProblemType(SVD,SVDProblemType);
381: SLEPC_EXTERN PetscErrorCode SVDGetProblemType(SVD,SVDProblemType*);
382: SLEPC_EXTERN PetscErrorCode SVDIsGeneralized(SVD,PetscBool*);
383: SLEPC_EXTERN PetscErrorCode SVDIsHyperbolic(SVD,PetscBool*);
384: SLEPC_EXTERN PetscErrorCode SVDSetOperators(SVD,Mat,Mat);
385: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDSetOperators()", ) static inline PetscErrorCode SVDSetOperator(SVD svd,Mat A) {return SVDSetOperators(svd,A,PETSC_NULLPTR);}
386: SLEPC_EXTERN PetscErrorCode SVDGetOperators(SVD,Mat*,Mat*);
387: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "SVDGetOperators()", ) static inline PetscErrorCode SVDGetOperator(SVD svd,Mat *A) {return SVDGetOperators(svd,A,PETSC_NULLPTR);}
388: SLEPC_EXTERN PetscErrorCode SVDSetSignature(SVD,Vec);
389: SLEPC_EXTERN PetscErrorCode SVDGetSignature(SVD,Vec);
390: SLEPC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec[],PetscInt,Vec[]);
391: 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);}
392: 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);}
393: SLEPC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
394: SLEPC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
395: SLEPC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
396: SLEPC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
397: SLEPC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
398: SLEPC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
399: SLEPC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
400: SLEPC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
401: SLEPC_EXTERN PetscErrorCode SVDSetThreshold(SVD,PetscReal,PetscBool);
402: SLEPC_EXTERN PetscErrorCode SVDGetThreshold(SVD,PetscReal*,PetscBool*);
403: SLEPC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
404: SLEPC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char[]);
405: SLEPC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char[]);
406: SLEPC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
407: SLEPC_EXTERN PetscErrorCode SVDSetDSType(SVD);
408: SLEPC_EXTERN PetscErrorCode SVDSetUp(SVD);
409: SLEPC_EXTERN PetscErrorCode SVDSolve(SVD);
410: SLEPC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
411: SLEPC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
412: SLEPC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
413: SLEPC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
414: SLEPC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
415: 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);}
416: 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);}
417: SLEPC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
418: SLEPC_EXTERN PetscErrorCode SVDViewFromOptions(SVD,PetscObject,const char[]);
419: SLEPC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
420: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "SVDErrorView()", ) static inline PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
421: SLEPC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
422: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonView(SVD,PetscViewer);
423: SLEPC_EXTERN PetscErrorCode SVDConvergedReasonViewFromOptions(SVD);
424: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonView()", ) static inline PetscErrorCode SVDReasonView(SVD svd,PetscViewer v) {return SVDConvergedReasonView(svd,v);}
425: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SVDConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SVDReasonViewFromOptions(SVD svd) {return SVDConvergedReasonViewFromOptions(svd);}
426: SLEPC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
427: SLEPC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
428: SLEPC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
429: SLEPC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
430: SLEPC_EXTERN PetscErrorCode SVDDestroy(SVD*);
431: SLEPC_EXTERN PetscErrorCode SVDReset(SVD);
432: SLEPC_EXTERN PetscErrorCode SVDSetWorkVecs(SVD,PetscInt,PetscInt);
433: SLEPC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
434: SLEPC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);

436: /*S
437:    SVDMonitorFn - A function prototype for functions provided to `SVDMonitorSet()`.

439:    Calling Sequence:
440: +  svd    - the singular value solver context
441: .  its    - iteration number
442: .  nconv  - number of converged singular triplets
443: .  sigma  - singular values
444: .  errest - relative error estimates for each singular triplet
445: .  nest   - number of error estimates
446: -  ctx    - optional monitoring context, as provided with `SVDMonitorSet()`

448:    Level: intermediate

450: .seealso: [](ch:svd), `SVDMonitorSet()`
451: S*/
452: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorFn(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,void *ctx);

454: /*S
455:    SVDMonitorRegisterFn - A function prototype for functions provided to `SVDMonitorRegister()`.

457:    Calling Sequence:
458: +  svd    - the singular value solver context
459: .  its    - iteration number
460: .  nconv  - number of converged singular triplets
461: .  sigma  - singular values
462: .  errest - relative error estimates for each singular triplet
463: .  nest   - number of error estimates
464: -  ctx    - `PetscViewerAndFormat` object

466:    Level: advanced

468:    Note:
469:    This is an `SVDMonitorFn` specialized for a context of `PetscViewerAndFormat`.

471: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorRegister()`, `SVDMonitorFn`, `SVDMonitorRegisterCreateFn`, `SVDMonitorRegisterDestroyFn`
472: S*/
473: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterFn(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *ctx);

475: /*S
476:    SVDMonitorRegisterCreateFn - A function prototype for functions that do
477:    the creation when provided to `SVDMonitorRegister()`.

479:    Calling Sequence:
480: +  viewer - the viewer to be used with the `SVDMonitorRegisterFn`
481: .  format - the format of the viewer
482: .  ctx    - a context for the monitor
483: -  result - a `PetscViewerAndFormat` object

485:    Level: advanced

487: .seealso: [](ch:svd), `SVDMonitorRegisterFn`, `SVDMonitorSet()`, `SVDMonitorRegister()`, `SVDMonitorFn`, `SVDMonitorRegisterDestroyFn`
488: S*/
489: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **result);

491: /*S
492:    SVDMonitorRegisterDestroyFn - A function prototype for functions that do the after
493:    use destruction when provided to `SVDMonitorRegister()`.

495:    Calling Sequence:
496: .  vf - a `PetscViewerAndFormat` object to be destroyed, including any context

498:    Level: advanced

500: .seealso: [](ch:svd), `SVDMonitorRegisterFn`, `SVDMonitorSet()`, `SVDMonitorRegister()`, `SVDMonitorFn`, `SVDMonitorRegisterCreateFn`
501: S*/
502: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDMonitorRegisterDestroyFn(PetscViewerAndFormat **result);

504: SLEPC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal[],PetscReal[],PetscInt);
505: SLEPC_EXTERN PetscErrorCode SVDMonitorSet(SVD,SVDMonitorFn,void*,PetscCtxDestroyFn*);
506: SLEPC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
507: SLEPC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void*);

509: SLEPC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char[],const char[],void*,PetscBool);
510: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorFirst;
511: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorFirstDrawLG;
512: SLEPC_EXTERN SVDMonitorRegisterCreateFn  SVDMonitorFirstDrawLGCreate;
513: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorAll;
514: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorAllDrawLG;
515: SLEPC_EXTERN SVDMonitorRegisterCreateFn  SVDMonitorAllDrawLGCreate;
516: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorConverged;
517: SLEPC_EXTERN SVDMonitorRegisterCreateFn  SVDMonitorConvergedCreate;
518: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorConvergedDrawLG;
519: SLEPC_EXTERN SVDMonitorRegisterCreateFn  SVDMonitorConvergedDrawLGCreate;
520: SLEPC_EXTERN SVDMonitorRegisterDestroyFn SVDMonitorConvergedDestroy;
521: SLEPC_EXTERN SVDMonitorRegisterFn        SVDMonitorConditioning;

523: SLEPC_EXTERN PetscFunctionList SVDList;
524: SLEPC_EXTERN PetscFunctionList SVDMonitorList;
525: SLEPC_EXTERN PetscFunctionList SVDMonitorCreateList;
526: SLEPC_EXTERN PetscFunctionList SVDMonitorDestroyList;
527: SLEPC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
528: SLEPC_EXTERN PetscErrorCode SVDMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,SVDMonitorRegisterFn*,SVDMonitorRegisterCreateFn*,SVDMonitorRegisterDestroyFn*);

530: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
531: SLEPC_EXTERN PetscErrorCode SVDReallocateSolution(SVD,PetscInt);

533: /*S
534:    SVDConvergenceTestFn - A prototype of an `SVD` convergence test function that
535:    would be passed to `SVDSetConvergenceTestFunction()`.

537:    Calling Sequence:
538: +  svd    - the singular value solver context
539: .  sigma  - computed singular value
540: .  res    - residual norm associated to the singular triplet
541: .  errest - [output] computed error estimate
542: -  ctx    - optional convergence context, as set by `SVDSetConvergenceTestFunction()`

544:    Level: advanced

546: .seealso: [](ch:svd), `SVDSetConvergenceTest()`, `SVDSetConvergenceTestFunction()`
547: S*/
548: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDConvergenceTestFn(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx);

550: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
551: SLEPC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
552: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedAbsolute;
553: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedRelative;
554: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedNorm;
555: SLEPC_EXTERN SVDConvergenceTestFn SVDConvergedMaxIt;
556: SLEPC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,SVDConvergenceTestFn*,void*,PetscCtxDestroyFn*);

558: /*S
559:    SVDStoppingTestFn - A prototype of an `SVD` stopping test function that would
560:    be passed to `SVDSetStoppingTestFunction()`.

562:    Calling Sequence:
563: +  svd    - the singular value solver context
564: .  its    - current number of iterations
565: .  max_it - maximum number of iterations
566: .  nconv  - number of currently converged singular triplets
567: .  nsv    - number of requested singular triplets
568: .  reason - [output] result of the stopping test
569: -  ctx    - optional stopping context, as set by `SVDSetStoppingTestFunction()`

571:    Note:
572:    A positive value of `reason` indicates that the iteration has finished successfully
573:    (converged), and a negative value indicates an error condition (diverged). If
574:    the iteration needs to be continued, `reason` must be set to `SVD_CONVERGED_ITERATING`
575:    (zero).

577:    Level: advanced

579: .seealso: [](ch:svd), `SVDSetStoppingTest()`, `SVDSetStoppingTestFunction()`
580: S*/
581: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SVDStoppingTestFn(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx);

583: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
584: SLEPC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
585: SLEPC_EXTERN SVDStoppingTestFn SVDStoppingBasic;
586: SLEPC_EXTERN SVDStoppingTestFn SVDStoppingThreshold;
587: SLEPC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,SVDStoppingTestFn*,void*,PetscCtxDestroyFn*);

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

591: SLEPC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
592: SLEPC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
593: SLEPC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
594: SLEPC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);

596: SLEPC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
597: SLEPC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
598: SLEPC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
599: SLEPC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);

601: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
602: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

604: /*E
605:    SVDTRLanczosGBidiag - The choice of bidiagonalization for the `SVDTRLANCZOS` GSVD solver.

607:    Values:
608: +  `SVD_TRLANCZOS_GBIDIAG_SINGLE` - single bidiagonalization ($Q_A$)
609: .  `SVD_TRLANCZOS_GBIDIAG_UPPER`  - joint bidiagonalization, both $Q_A$ and $Q_B$ in upper bidiagonal form
610: -  `SVD_TRLANCZOS_GBIDIAG_LOWER`  - joint bidiagonalization, $Q_A$ lower bidiagonal, $Q_B$ upper bidiagonal

612:    Note:
613:    The different variants are described in {cite:p}`Alv24`.

615:    Level: advanced

617: .seealso: [](ch:svd), `SVDTRLanczosSetGBidiag()`, `SVDTRLanczosGetGBidiag()`
618: E*/
619: typedef enum {
620:   SVD_TRLANCZOS_GBIDIAG_SINGLE,
621:   SVD_TRLANCZOS_GBIDIAG_UPPER,
622:   SVD_TRLANCZOS_GBIDIAG_LOWER
623: } SVDTRLanczosGBidiag;
624: SLEPC_EXTERN const char *SVDTRLanczosGBidiags[];

626: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetGBidiag(SVD,SVDTRLanczosGBidiag);
627: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetGBidiag(SVD,SVDTRLanczosGBidiag*);
628: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
629: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
630: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetKSP(SVD,KSP);
631: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetKSP(SVD,KSP*);
632: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetRestart(SVD,PetscReal);
633: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetRestart(SVD,PetscReal*);
634: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetLocking(SVD,PetscBool);
635: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetLocking(SVD,PetscBool*);
636: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD,PetscBool);
637: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD,PetscBool*);
638: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetScale(SVD,PetscReal);
639: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetScale(SVD,PetscReal*);

641: /*E
642:    SVDPRIMMEMethod - The SVD method selected in the PRIMME library.

644:    Note:
645:    See the documentation of PRIMME {cite:p}`Sta10` for a description of the methods.

647:    Level: advanced

649: .seealso: [](ch:svd), `SVDPRIMMESetMethod()`, `SVDPRIMMEGetMethod()`
650: E*/
651: typedef enum { SVD_PRIMME_HYBRID          = 1,
652:                SVD_PRIMME_NORMALEQUATIONS = 2,
653:                SVD_PRIMME_AUGMENTED       = 3 } SVDPRIMMEMethod;
654: SLEPC_EXTERN const char *SVDPRIMMEMethods[];

656: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
657: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
658: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
659: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);

661: /*E
662:    SVDKSVDEigenMethod - The method to solve the eigenproblem within the KSVD library.

664:    Note:
665:    See the documentation of KSVD {cite:p}`Suk19` for a description of the methods.

667:    Level: advanced

669: .seealso: [](ch:svd), `SVDKSVDSetEigenMethod()`, `SVDKSVDGetEigenMethod()`
670: E*/
671: typedef enum { SVD_KSVD_EIGEN_MRRR = 1,
672:                SVD_KSVD_EIGEN_DC   = 2,
673:                SVD_KSVD_EIGEN_ELPA = 3 } SVDKSVDEigenMethod;
674: SLEPC_EXTERN const char *SVDKSVDEigenMethods[];

676: /*E
677:    SVDKSVDPolarMethod - The method to compute the polar decomposition within the KSVD library.

679:    Note:
680:    See the documentation of KSVD {cite:p}`Suk19` for a description of the methods.

682:    Level: advanced

684: .seealso: [](ch:svd), `SVDKSVDSetPolarMethod()`, `SVDKSVDGetPolarMethod()`
685: E*/
686: typedef enum { SVD_KSVD_POLAR_QDWH   = 1,
687:                SVD_KSVD_POLAR_ZOLOPD = 2 } SVDKSVDPolarMethod;
688: SLEPC_EXTERN const char *SVDKSVDPolarMethods[];

690: SLEPC_EXTERN PetscErrorCode SVDKSVDSetEigenMethod(SVD,SVDKSVDEigenMethod);
691: SLEPC_EXTERN PetscErrorCode SVDKSVDGetEigenMethod(SVD,SVDKSVDEigenMethod*);
692: SLEPC_EXTERN PetscErrorCode SVDKSVDSetPolarMethod(SVD,SVDKSVDPolarMethod);
693: SLEPC_EXTERN PetscErrorCode SVDKSVDGetPolarMethod(SVD,SVDKSVDPolarMethod*);