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