Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1251 slepc 1
/*
2
   User interface for the SLEPC singular value solvers.
3
*/
4
#if !defined(__SLEPCSVD_H)
5
#define __SLEPCSVD_H
6
#include "slepc.h"
1253 slepc 7
#include "slepceps.h"
1251 slepc 8
PETSC_EXTERN_CXX_BEGIN
9
 
10
extern PetscCookie SVD_COOKIE;
11
 
12
/*S
13
     SVD - Abstract SLEPc object that manages all the singular value
14
     problem solvers.
15
 
16
   Level: beginner
17
 
18
.seealso:  SVDCreate()
19
S*/
20
typedef struct _p_SVD* SVD;
21
 
1364 slepc 22
/*E
23
    SVDType - String with the name of a SLEPc singular value solver
24
 
25
   Level: beginner
26
 
27
.seealso: SVDSetType(), SVD
28
E*/
1251 slepc 29
#define SVDType const char*
1324 slepc 30
#define SVDCROSS       "cross"
31
#define SVDCYCLIC      "cyclic"
1251 slepc 32
#define SVDLAPACK      "lapack"
1268 slepc 33
#define SVDLANCZOS     "lanczos"
1298 slepc 34
#define SVDTRLANCZOS   "trlanczos"
1251 slepc 35
 
1364 slepc 36
/*E
37
    SVDTransposeMode - determines how to handle the transpose of the matrix
38
 
39
    Level: advanced
40
 
41
.seealso: SVDSetTransposeMode(), SVDGetTransposeMode()
42
E*/
1283 slepc 43
typedef enum { SVD_TRANSPOSE_EXPLICIT, SVD_TRANSPOSE_IMPLICIT } SVDTransposeMode;
1255 slepc 44
 
1364 slepc 45
/*E
46
    SVDWhich - determines whether largest or smallest singular triplets
47
    are to be computed
48
 
49
    Level: intermediate
50
 
51
.seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
52
E*/
1272 slepc 53
typedef enum { SVD_LARGEST, SVD_SMALLEST } SVDWhich;
54
 
1364 slepc 55
/*E
56
    SVDConvergedReason - reason a singular value solver was said to
57
         have converged or diverged
58
 
59
   Level: beginner
60
 
61
.seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
62
E*/
1283 slepc 63
typedef enum {/* converged */
64
              SVD_CONVERGED_TOL                =  2,
65
              /* diverged */
66
              SVD_DIVERGED_ITS                 = -3,
67
              SVD_DIVERGED_BREAKDOWN           = -4,
68
              SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;
69
 
1251 slepc 70
EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
1345 slepc 71
EXTERN PetscErrorCode SVDSetIP(SVD,IP);
72
EXTERN PetscErrorCode SVDGetIP(SVD,IP*);
1251 slepc 73
EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
74
EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
75
EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
1268 slepc 76
EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
1272 slepc 77
EXTERN PetscErrorCode SVDSetInitialVector(SVD,Vec);
78
EXTERN PetscErrorCode SVDGetInitialVector(SVD,Vec*);
1268 slepc 79
EXTERN PetscErrorCode SVDSetTransposeMode(SVD,SVDTransposeMode);
80
EXTERN PetscErrorCode SVDGetTransposeMode(SVD,SVDTransposeMode*);
1272 slepc 81
EXTERN PetscErrorCode SVDSetDimensions(SVD,int,int);
82
EXTERN PetscErrorCode SVDGetDimensions(SVD,int*,int*);
83
EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,int);
84
EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,int*);
85
EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
86
EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
1251 slepc 87
EXTERN PetscErrorCode SVDSetFromOptions(SVD);
1309 slepc 88
EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
89
EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
90
EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
1251 slepc 91
EXTERN PetscErrorCode SVDSetUp(SVD);
92
EXTERN PetscErrorCode SVDSolve(SVD);
1283 slepc 93
EXTERN PetscErrorCode SVDGetIterationNumber(SVD,int*);
1337 slepc 94
EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
1251 slepc 95
EXTERN PetscErrorCode SVDGetConverged(SVD,int*);
96
EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,int,PetscReal*,Vec,Vec);
1257 slepc 97
EXTERN PetscErrorCode SVDComputeResidualNorms(SVD,int,PetscReal*,PetscReal*);
1317 slepc 98
EXTERN PetscErrorCode SVDComputeRelativeError(SVD,int,PetscReal*);
1305 slepc 99
EXTERN PetscErrorCode SVDGetOperationCounters(SVD,int*,int*);
1251 slepc 100
EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
101
EXTERN PetscErrorCode SVDDestroy(SVD);
102
EXTERN PetscErrorCode SVDInitializePackage(char*);
103
 
1331 slepc 104
EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,int,int,PetscReal*,PetscReal*,int,void*),
1288 slepc 105
                                    void*,PetscErrorCode (*monitordestroy)(void*));
1331 slepc 106
EXTERN PetscErrorCode SVDMonitorCancel(SVD);
1288 slepc 107
EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
1331 slepc 108
EXTERN PetscErrorCode SVDMonitorDefault(SVD,int,int,PetscReal*,PetscReal*,int,void*);
109
EXTERN PetscErrorCode SVDMonitorLG(SVD,int,int,PetscReal*,PetscReal*,int,void*);
1288 slepc 110
 
1298 slepc 111
EXTERN PetscErrorCode SVDDense(int,int,PetscScalar*,PetscReal*,PetscScalar*,PetscScalar*);
112
 
1324 slepc 113
EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
114
EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
115
 
116
EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscTruth);
117
EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscTruth*);
118
EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
119
EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
120
 
1359 slepc 121
EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscTruth);
1298 slepc 122
 
1359 slepc 123
EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscTruth);
1298 slepc 124
 
1337 slepc 125
EXTERN PetscErrorCode SVDDense(int,int,PetscScalar*,PetscReal*,PetscScalar*,PetscScalar*);
126
 
1251 slepc 127
#endif