Actual source code: svdimpl.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc. See the README file for conditions of use
7: and additional information.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: #ifndef _SVDIMPL
12: #define _SVDIMPL
14: #include slepcsvd.h
15: #include slepcip.h
20: typedef struct _SVDOps *SVDOps;
22: struct _SVDOps {
23: int (*solve)(SVD);
24: int (*setup)(SVD);
25: int (*setfromoptions)(SVD);
26: int (*publishoptions)(SVD);
27: int (*destroy)(SVD);
28: int (*view)(SVD,PetscViewer);
29: };
31: /*
32: Maximum number of monitors you can run with a single SVD
33: */
34: #define MAXSVDMONITORS 5
36: /*
37: Defines the SVD data structure.
38: */
39: struct _p_SVD {
40: PETSCHEADER(struct _SVDOps);
41: Mat OP; /* problem matrix */
42: Mat A; /* problem matrix (m>n) */
43: Mat AT; /* transposed matrix */
44: SVDTransposeMode transmode; /* transpose mode */
45: PetscReal *sigma; /* singular values */
46: Vec *U,*V; /* left and right singular vectors */
47: Vec vec_initial; /* initial vector */
48: int n; /* maximun size of descomposition */
49: SVDWhich which; /* which singular values are computed */
50: int nconv; /* number of converged values */
51: int nsv; /* number of requested values */
52: int ncv; /* basis size */
53: int its; /* iteration counter */
54: int max_it; /* max iterations */
55: PetscReal tol; /* tolerance */
56: PetscReal *errest; /* error estimates */
57: void *data; /* placeholder for misc stuff associated
58: with a particular solver */
59: int setupcalled;
60: SVDConvergedReason reason;
61: IP ip;
62:
63: int (*monitor[MAXSVDMONITORS])(SVD,int,int,PetscReal*,PetscReal*,int,void*);
64: int (*monitordestroy[MAXSVDMONITORS])(void*);
65: void *monitorcontext[MAXSVDMONITORS];
66: int numbermonitors;
67:
68: int matvecs;
69: };
71: EXTERN PetscErrorCode SVDRegisterAll(char *);
73: #define SVDMonitor(svd,it,nconv,sigma,errest,nest) \
74: { int _ierr,_i,_im = svd->numbermonitors; \
75: for ( _i=0; _i<_im; _i++ ) {\
76: _ierr=(*svd->monitor[_i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[_i]);\
77: CHKERRQ(_ierr); \
78: } \
79: }
81: #endif
82: EXTERN PetscErrorCode SVDDestroy_Default(SVD);
83: EXTERN PetscErrorCode SVDMatMult(SVD,PetscTruth,Vec,Vec);
84: EXTERN PetscErrorCode SVDMatGetVecs(SVD,Vec*,Vec*);
85: EXTERN PetscErrorCode SVDMatGetSize(SVD,PetscInt*,PetscInt*);
86: EXTERN PetscErrorCode SVDMatGetLocalSize(SVD,PetscInt*,PetscInt*);
87: EXTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,Vec*,Vec,Vec*,int,int,PetscScalar*,Vec,Vec);