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.
1376 slepc 3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
6
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
1376 slepc 7
 
1672 slepc 8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1251 slepc 22
*/
1376 slepc 23
 
1251 slepc 24
#if !defined(__SLEPCSVD_H)
25
#define __SLEPCSVD_H
26
#include "slepc.h"
1253 slepc 27
#include "slepceps.h"
1251 slepc 28
PETSC_EXTERN_CXX_BEGIN
29
 
30
extern PetscCookie SVD_COOKIE;
31
 
32
/*S
33
     SVD - Abstract SLEPc object that manages all the singular value
34
     problem solvers.
35
 
36
   Level: beginner
37
 
38
.seealso:  SVDCreate()
39
S*/
40
typedef struct _p_SVD* SVD;
41
 
1364 slepc 42
/*E
43
    SVDType - String with the name of a SLEPc singular value solver
44
 
45
   Level: beginner
46
 
47
.seealso: SVDSetType(), SVD
48
E*/
1502 slepc 49
#define SVDType        char*
1324 slepc 50
#define SVDCROSS       "cross"
51
#define SVDCYCLIC      "cyclic"
1251 slepc 52
#define SVDLAPACK      "lapack"
1268 slepc 53
#define SVDLANCZOS     "lanczos"
1298 slepc 54
#define SVDTRLANCZOS   "trlanczos"
1251 slepc 55
 
1364 slepc 56
/*E
57
    SVDTransposeMode - determines how to handle the transpose of the matrix
58
 
59
    Level: advanced
60
 
61
.seealso: SVDSetTransposeMode(), SVDGetTransposeMode()
62
E*/
1283 slepc 63
typedef enum { SVD_TRANSPOSE_EXPLICIT, SVD_TRANSPOSE_IMPLICIT } SVDTransposeMode;
1255 slepc 64
 
1364 slepc 65
/*E
66
    SVDWhich - determines whether largest or smallest singular triplets
67
    are to be computed
68
 
69
    Level: intermediate
70
 
71
.seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
72
E*/
1272 slepc 73
typedef enum { SVD_LARGEST, SVD_SMALLEST } SVDWhich;
74
 
1364 slepc 75
/*E
76
    SVDConvergedReason - reason a singular value solver was said to
77
         have converged or diverged
78
 
79
   Level: beginner
80
 
81
.seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
82
E*/
1283 slepc 83
typedef enum {/* converged */
84
              SVD_CONVERGED_TOL                =  2,
85
              /* diverged */
86
              SVD_DIVERGED_ITS                 = -3,
87
              SVD_DIVERGED_BREAKDOWN           = -4,
88
              SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;
89
 
1251 slepc 90
EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
1345 slepc 91
EXTERN PetscErrorCode SVDSetIP(SVD,IP);
92
EXTERN PetscErrorCode SVDGetIP(SVD,IP*);
1502 slepc 93
EXTERN PetscErrorCode SVDSetType(SVD,const SVDType);
1501 slepc 94
EXTERN PetscErrorCode SVDGetType(SVD,const SVDType*);
1251 slepc 95
EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
1268 slepc 96
EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
1272 slepc 97
EXTERN PetscErrorCode SVDSetInitialVector(SVD,Vec);
98
EXTERN PetscErrorCode SVDGetInitialVector(SVD,Vec*);
1268 slepc 99
EXTERN PetscErrorCode SVDSetTransposeMode(SVD,SVDTransposeMode);
100
EXTERN PetscErrorCode SVDGetTransposeMode(SVD,SVDTransposeMode*);
1575 slepc 101
EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
102
EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
1504 slepc 103
EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
104
EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
1272 slepc 105
EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
106
EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
1251 slepc 107
EXTERN PetscErrorCode SVDSetFromOptions(SVD);
1309 slepc 108
EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
109
EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
110
EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
1251 slepc 111
EXTERN PetscErrorCode SVDSetUp(SVD);
112
EXTERN PetscErrorCode SVDSolve(SVD);
1504 slepc 113
EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
1337 slepc 114
EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
1504 slepc 115
EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
116
EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
117
EXTERN PetscErrorCode SVDComputeResidualNorms(SVD,PetscInt,PetscReal*,PetscReal*);
118
EXTERN PetscErrorCode SVDComputeRelativeError(SVD,PetscInt,PetscReal*);
119
EXTERN PetscErrorCode SVDGetOperationCounters(SVD,PetscInt*,PetscInt*);
1251 slepc 120
EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
121
EXTERN PetscErrorCode SVDDestroy(SVD);
122
EXTERN PetscErrorCode SVDInitializePackage(char*);
123
 
1504 slepc 124
EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),
1288 slepc 125
                                    void*,PetscErrorCode (*monitordestroy)(void*));
1331 slepc 126
EXTERN PetscErrorCode SVDMonitorCancel(SVD);
1288 slepc 127
EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
1504 slepc 128
EXTERN PetscErrorCode SVDMonitorDefault(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
129
EXTERN PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
1288 slepc 130
 
1504 slepc 131
EXTERN PetscErrorCode SVDDense(PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscScalar*,PetscScalar*);
1298 slepc 132
 
1324 slepc 133
EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
134
EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
135
 
136
EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscTruth);
137
EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscTruth*);
138
EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
139
EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
140
 
1359 slepc 141
EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscTruth);
1298 slepc 142
 
1359 slepc 143
EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscTruth);
1298 slepc 144
 
1504 slepc 145
EXTERN PetscErrorCode SVDRegister(const char*,const char*,const char*,PetscErrorCode(*)(SVD));
1389 slepc 146
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
147
#define SVDRegisterDynamic(a,b,c,d) SVDRegister(a,b,c,0)
148
#else
149
#define SVDRegisterDynamic(a,b,c,d) SVDRegister(a,b,c,d)
1251 slepc 150
#endif
1389 slepc 151
EXTERN PetscErrorCode SVDRegisterDestroy(void);
152
 
1410 slepc 153
PETSC_EXTERN_CXX_END
1389 slepc 154
#endif