Actual source code: svdimpl.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
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: #if !defined(_SVDIMPL)
23: #define _SVDIMPL
25: #include <slepcsvd.h>
26: #include <slepc-private/slepcimpl.h>
28: PETSC_EXTERN PetscLogEvent SVD_SetUp,SVD_Solve;
30: typedef struct _SVDOps *SVDOps;
32: struct _SVDOps {
33: PetscErrorCode (*solve)(SVD);
34: PetscErrorCode (*setup)(SVD);
35: PetscErrorCode (*setfromoptions)(SVD);
36: PetscErrorCode (*publishoptions)(SVD);
37: PetscErrorCode (*destroy)(SVD);
38: PetscErrorCode (*reset)(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 *IS,*ISL; /* placeholder for references to user-provided initial space */
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 nini,ninil; /* number of initial vectors (negative means not copied yet) */
67: PetscInt its; /* iteration counter */
68: PetscInt max_it; /* max iterations */
69: PetscReal tol; /* tolerance */
70: PetscReal *errest; /* error estimates */
71: PetscRandom rand; /* random number generator */
72: Vec tl,tr; /* template vectors */
73: void *data; /* placeholder for misc stuff associated
74: with a particular solver */
75: PetscInt setupcalled;
76: SVDConvergedReason reason;
77: IP ip; /* innerproduct object */
78: DS ds; /* direct solver object */
79: PetscBool trackall;
80: PetscInt matvecs;
82: PetscErrorCode (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
83: PetscErrorCode (*monitordestroy[MAXSVDMONITORS])(void**);
84: void *monitorcontext[MAXSVDMONITORS];
85: PetscInt numbermonitors;
86: };
88: PETSC_INTERN PetscErrorCode SVDMatMult(SVD,PetscBool,Vec,Vec);
89: PETSC_INTERN PetscErrorCode SVDMatGetVecs(SVD,Vec*,Vec*);
90: PETSC_INTERN PetscErrorCode SVDMatGetSize(SVD,PetscInt*,PetscInt*);
91: PETSC_INTERN PetscErrorCode SVDMatGetLocalSize(SVD,PetscInt*,PetscInt*);
92: PETSC_INTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,Vec*,Vec,Vec*,PetscInt,PetscInt,PetscScalar*);
94: #endif