Actual source code: slepcsvd.h

  1: /*
  2:    User interface for the SLEPC singular value solvers. 

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  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.

 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.

 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/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 26:  #include slepcsys.h
 27:  #include slepceps.h

 29: PETSC_EXTERN PetscErrorCode SVDInitializePackage(const char[]);

 31: /*S
 32:      SVD - Abstract SLEPc object that manages all the singular value 
 33:      problem solvers.

 35:    Level: beginner

 37: .seealso:  SVDCreate()
 38: S*/
 39: typedef struct _p_SVD* SVD;

 41: /*E
 42:     SVDType - String with the name of a SLEPc singular value solver

 44:    Level: beginner

 46: .seealso: SVDSetType(), SVD
 47: E*/
 48: #define SVDType        char*
 49: #define SVDCROSS       "cross"
 50: #define SVDCYCLIC      "cyclic"
 51: #define SVDLAPACK      "lapack"
 52: #define SVDLANCZOS     "lanczos"
 53: #define SVDTRLANCZOS   "trlanczos"

 55: /* Logging support */
 56: PETSC_EXTERN PetscClassId SVD_CLASSID;

 58: /*E
 59:     SVDTransposeMode - determines how to handle the transpose of the matrix

 61:     Level: advanced

 63: .seealso: SVDSetTransposeMode(), SVDGetTransposeMode()
 64: E*/
 65: typedef enum { SVD_TRANSPOSE_EXPLICIT,
 66:                SVD_TRANSPOSE_IMPLICIT } SVDTransposeMode;

 68: /*E
 69:     SVDWhich - determines whether largest or smallest singular triplets
 70:     are to be computed

 72:     Level: intermediate

 74: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
 75: E*/
 76: typedef enum { SVD_LARGEST,
 77:                SVD_SMALLEST } SVDWhich;

 79: /*E
 80:     SVDConvergedReason - reason a singular value solver was said to 
 81:          have converged or diverged

 83:    Level: beginner

 85: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
 86: E*/
 87: typedef enum {/* converged */
 88:               SVD_CONVERGED_TOL                =  2,
 89:               /* diverged */
 90:               SVD_DIVERGED_ITS                 = -3,
 91:               SVD_DIVERGED_BREAKDOWN           = -4,
 92:               SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;

 94: PETSC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
 95: PETSC_EXTERN PetscErrorCode SVDSetIP(SVD,IP);
 96: PETSC_EXTERN PetscErrorCode SVDGetIP(SVD,IP*);
 97: PETSC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
 98: PETSC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
 99: PETSC_EXTERN PetscErrorCode SVDSetType(SVD,const SVDType);
100: PETSC_EXTERN PetscErrorCode SVDGetType(SVD,const SVDType*);
101: PETSC_EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
102: PETSC_EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
103: PETSC_EXTERN PetscErrorCode SVDSetInitialSpace(SVD,PetscInt,Vec*);
104: PETSC_EXTERN PetscErrorCode SVDSetTransposeMode(SVD,SVDTransposeMode);
105: PETSC_EXTERN PetscErrorCode SVDGetTransposeMode(SVD,SVDTransposeMode*);
106: PETSC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
107: PETSC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
108: PETSC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
109: PETSC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
110: PETSC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
111: PETSC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
112: PETSC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
113: PETSC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
114: PETSC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
115: PETSC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
116: PETSC_EXTERN PetscErrorCode SVDSetUp(SVD);
117: PETSC_EXTERN PetscErrorCode SVDSolve(SVD);
118: PETSC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
119: PETSC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
120: PETSC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
121: PETSC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
122: PETSC_EXTERN PetscErrorCode SVDComputeResidualNorms(SVD,PetscInt,PetscReal*,PetscReal*);
123: PETSC_EXTERN PetscErrorCode SVDComputeRelativeError(SVD,PetscInt,PetscReal*);
124: PETSC_EXTERN PetscErrorCode SVDGetOperationCounters(SVD,PetscInt*,PetscInt*);
125: PETSC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
126: PETSC_EXTERN PetscErrorCode SVDPrintSolution(SVD,PetscViewer);
127: PETSC_EXTERN PetscErrorCode SVDDestroy(SVD*);
128: PETSC_EXTERN PetscErrorCode SVDReset(SVD);

130: PETSC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
131: PETSC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
132: PETSC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
133: PETSC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
134: PETSC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
135: PETSC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
136: PETSC_EXTERN PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
137: PETSC_EXTERN PetscErrorCode SVDMonitorLGAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);

139: PETSC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
140: PETSC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);

142: PETSC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
143: PETSC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);

145: PETSC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
146: PETSC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
147: PETSC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
148: PETSC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);

150: PETSC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
151: PETSC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

153: PETSC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
154: PETSC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);

156: PETSC_EXTERN PetscFList SVDList;
157: PETSC_EXTERN PetscBool  SVDRegisterAllCalled;
158: PETSC_EXTERN PetscErrorCode SVDRegisterAll(const char[]);
159: PETSC_EXTERN PetscErrorCode SVDRegisterDestroy(void);
160: PETSC_EXTERN PetscErrorCode SVDRegister(const char[],const char[],const char[],PetscErrorCode(*)(SVD));

162: /*MC
163:    SVDRegisterDynamic - Adds a method to the singular value solver package.

165:    Synopsis:
166:    PetscErrorCode SVDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(SVD))

168:    Not Collective

170:    Input Parameters:
171: +  name_solver - name of a new user-defined solver
172: .  path - path (either absolute or relative) the library containing this solver
173: .  name_create - name of routine to create the solver context
174: -  routine_create - routine to create the solver context

176:    Notes:
177:    SVDRegisterDynamic() may be called multiple times to add several user-defined solvers.

179:    If dynamic libraries are used, then the fourth input argument (routine_create)
180:    is ignored.

182:    Sample usage:
183: .vb
184:    SVDRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
185:                "MySolverCreate",MySolverCreate);
186: .ve

188:    Then, your solver can be chosen with the procedural interface via
189: $     SVDSetType(svd,"my_solver")
190:    or at runtime via the option
191: $     -svd_type my_solver

193:    Level: advanced

195: .seealso: SVDRegisterDestroy(), SVDRegisterAll()
196: M*/
197: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
198: #define SVDRegisterDynamic(a,b,c,d) SVDRegister(a,b,c,0)
199: #else
200: #define SVDRegisterDynamic(a,b,c,d) SVDRegister(a,b,c,d)
201: #endif

203: #endif