Actual source code: svdimpl.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat 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>
28: extern PetscFList SVDList;
29: extern PetscLogEvent SVD_SetUp, SVD_Solve, SVD_Dense;
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 (*reset)(SVD);
40: PetscErrorCode (*view)(SVD,PetscViewer);
41: };
43: /*
44: Maximum number of monitors you can run with a single SVD
45: */
46: #define MAXSVDMONITORS 5
48: /*
49: Defines the SVD data structure.
50: */
51: struct _p_SVD {
52: PETSCHEADER(struct _SVDOps);
53: Mat OP; /* problem matrix */
54: Mat A; /* problem matrix (m>n) */
55: Mat AT; /* transposed matrix */
56: SVDTransposeMode transmode; /* transpose mode */
57: PetscReal *sigma; /* singular values */
58: PetscInt *perm; /* permutation for singular value ordering */
59: Vec *U,*V; /* left and right singular vectors */
60: Vec *IS; /* placeholder for references to user-provided initial space */
61: PetscInt n; /* maximun size of descomposition */
62: SVDWhich which; /* which singular values are computed */
63: PetscInt nconv; /* number of converged values */
64: PetscInt nsv; /* number of requested values */
65: PetscInt ncv; /* basis size */
66: PetscInt mpd; /* maximum dimension of projected problem */
67: PetscInt nini; /* number of initial vectors (negative means not copied yet) */
68: PetscInt its; /* iteration counter */
69: PetscInt max_it; /* max iterations */
70: PetscReal tol; /* tolerance */
71: PetscReal *errest; /* error estimates */
72: PetscRandom rand; /* random number generator */
73: Vec tl,tr; /* template vectors */
74: void *data; /* placeholder for misc stuff associated
75: with a particular solver */
76: PetscInt setupcalled;
77: SVDConvergedReason reason;
78: IP ip;
79: PetscBool trackall;
80:
81: PetscErrorCode (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
82: PetscErrorCode (*monitordestroy[MAXSVDMONITORS])(void**);
83: void *monitorcontext[MAXSVDMONITORS];
84: PetscInt numbermonitors;
85:
86: PetscInt matvecs;
87: };
89: extern PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
91: extern PetscErrorCode SVDMatMult(SVD,PetscBool,Vec,Vec);
92: extern PetscErrorCode SVDMatGetVecs(SVD,Vec*,Vec*);
93: extern PetscErrorCode SVDMatGetSize(SVD,PetscInt*,PetscInt*);
94: extern PetscErrorCode SVDMatGetLocalSize(SVD,PetscInt*,PetscInt*);
95: extern PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,Vec*,Vec,Vec*,PetscInt,PetscInt,PetscScalar*);
97: #endif