Actual source code: slepcsvd.h

slepc-3.15.0 2021-03-31
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_MAXIT,
 96:                SVD_CONV_USER } SVDConv;

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

101:     Level: advanced

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

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

112:    Level: intermediate

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

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

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

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

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

219: SLEPC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);

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

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

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

233: SLEPC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
234: SLEPC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

236: SLEPC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
237: SLEPC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);

239: /*E
240:     SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library

242:     Level: advanced

244: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
245: E*/
246: typedef enum { SVD_PRIMME_HYBRID=1,
247:                SVD_PRIMME_NORMALEQUATIONS,
248:                SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
249: SLEPC_EXTERN const char *SVDPRIMMEMethods[];

251: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
252: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
253: SLEPC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
254: SLEPC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);

256: #endif