Actual source code: svdimpl.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #ifndef _SVDIMPL
23: #define _SVDIMPL
25: #include slepcsvd.h
26: #include slepcip.h
31: typedef struct _SVDOps *SVDOps;
33: struct _SVDOps {
34: PetscErrorCode (*solve)(SVD);
35: PetscErrorCode (*setup)(SVD);
36: PetscErrorCode (*setfromoptions)(SVD);
37: PetscErrorCode (*publishoptions)(SVD);
38: PetscErrorCode (*destroy)(SVD);
39: PetscErrorCode (*view)(SVD,PetscViewer);
40: };
42: /*
43: Maximum number of monitors you can run with a single SVD
44: */
45: #define MAXSVDMONITORS 5
47: /*
48: Defines the SVD data structure.
49: */
50: struct _p_SVD {
51: PETSCHEADER(struct _SVDOps);
52: Mat OP; /* problem matrix */
53: Mat A; /* problem matrix (m>n) */
54: Mat AT; /* transposed matrix */
55: SVDTransposeMode transmode; /* transpose mode */
56: PetscReal *sigma; /* singular values */
57: PetscInt *perm; /* permutation for singular value ordering */
58: Vec *U,*V; /* left and right singular vectors */
59: Vec vec_initial; /* initial vector */
60: PetscInt n; /* maximun size of descomposition */
61: SVDWhich which; /* which singular values are computed */
62: PetscInt nconv; /* number of converged values */
63: PetscInt nsv; /* number of requested values */
64: PetscInt ncv; /* basis size */
65: PetscInt mpd; /* maximum dimension of projected problem */
66: PetscInt its; /* iteration counter */
67: PetscInt max_it; /* max iterations */
68: PetscReal tol; /* tolerance */
69: PetscReal *errest; /* error estimates */
70: void *data; /* placeholder for misc stuff associated
71: with a particular solver */
72: PetscInt setupcalled;
73: SVDConvergedReason reason;
74: IP ip;
75:
76: PetscErrorCode (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
77: PetscErrorCode (*monitordestroy[MAXSVDMONITORS])(void*);
78: void *monitorcontext[MAXSVDMONITORS];
79: PetscInt numbermonitors;
80:
81: PetscInt matvecs;
82: };
84: EXTERN PetscErrorCode SVDRegisterAll(char *);
86: #define SVDMonitor(svd,it,nconv,sigma,errest,nest) \
87: { PetscErrorCode _ierr; PetscInt _i,_im = svd->numbermonitors; \
88: for ( _i=0; _i<_im; _i++ ) {\
89: _ierr=(*svd->monitor[_i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[_i]);\
90: CHKERRQ(_ierr); \
91: } \
92: }
94: #endif
95: EXTERN PetscErrorCode SVDDestroy_Default(SVD);
96: EXTERN PetscErrorCode SVDMatMult(SVD,PetscTruth,Vec,Vec);
97: EXTERN PetscErrorCode SVDMatGetVecs(SVD,Vec*,Vec*);
98: EXTERN PetscErrorCode SVDMatGetSize(SVD,PetscInt*,PetscInt*);
99: EXTERN PetscErrorCode SVDMatGetLocalSize(SVD,PetscInt*,PetscInt*);
100: EXTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,Vec*,Vec,Vec*,PetscInt,PetscInt,PetscScalar*,Vec,Vec);