Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1521 slepc 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 3
   SLEPc - Scalable Library for Eigenvalue Problem Computations
4
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
1521 slepc 5
 
1672 slepc 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.
11
 
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.
16
 
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/>.
1521 slepc 19
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20
*/
21
 
22
#ifndef _SVDIMPL
23
#define _SVDIMPL
24
 
25
#include "slepcsvd.h"
26
#include "slepcip.h"
27
 
28
extern PetscFList SVDList;
29
extern PetscLogEvent SVD_SetUp, SVD_Solve, SVD_Dense;
30
 
31
typedef struct _SVDOps *SVDOps;
32
 
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
};
41
 
42
/*
43
     Maximum number of monitors you can run with a single SVD
44
*/
45
#define MAXSVDMONITORS 5
46
 
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 */
1603 slepc 57
  PetscInt         *perm;       /* permutation for singular value ordering */
1521 slepc 58
  Vec              *U,*V;       /* left and right singular vectors */
1952 jroman 59
  Vec              *IS;         /* placeholder for references to user-provided initial space */
1521 slepc 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 */
1575 slepc 65
  PetscInt         mpd;         /* maximum dimension of projected problem */
1952 jroman 66
  PetscInt         nini;        /* number of initial vectors (negative means not copied yet) */
1521 slepc 67
  PetscInt         its;         /* iteration counter */
68
  PetscInt         max_it;      /* max iterations */
69
  PetscReal        tol;         /* tolerance */
70
  PetscReal        *errest;     /* error estimates */
2027 jroman 71
  PetscRandom      rand;        /* random number generator */
1521 slepc 72
  void             *data;       /* placeholder for misc stuff associated
73
                                   with a particular solver */
74
  PetscInt         setupcalled;
75
  SVDConvergedReason reason;
76
  IP               ip;
77
 
78
  PetscErrorCode  (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
79
  PetscErrorCode  (*monitordestroy[MAXSVDMONITORS])(void*);
80
  void            *monitorcontext[MAXSVDMONITORS];
81
  PetscInt        numbermonitors;
82
 
83
  PetscInt        matvecs;
84
};
85
 
86
EXTERN PetscErrorCode SVDRegisterAll(char *);
1853 antodo 87
EXTERN PetscErrorCode SVDInitializePackage(char*);
88
EXTERN PetscErrorCode SVDFinalizePackage(void);
1521 slepc 89
 
90
#define SVDMonitor(svd,it,nconv,sigma,errest,nest) \
91
        { PetscErrorCode _ierr; PetscInt _i,_im = svd->numbermonitors; \
92
          for ( _i=0; _i<_im; _i++ ) {\
93
            _ierr=(*svd->monitor[_i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[_i]);\
94
            CHKERRQ(_ierr); \
95
          } \
96
        }
97
 
1951 jroman 98
/* context for SVDMonitorConverged */
99
typedef struct {
100
  PetscViewerASCIIMonitor viewer;
101
  PetscInt oldnconv;
102
} SVDMONITOR_CONV;
103
EXTERN PetscErrorCode SVDMonitorDestroy_Converged(SVDMONITOR_CONV*);
104
 
1521 slepc 105
EXTERN PetscErrorCode SVDDestroy_Default(SVD);
106
EXTERN PetscErrorCode SVDMatMult(SVD,PetscTruth,Vec,Vec);
107
EXTERN PetscErrorCode SVDMatGetVecs(SVD,Vec*,Vec*);
108
EXTERN PetscErrorCode SVDMatGetSize(SVD,PetscInt*,PetscInt*);
109
EXTERN PetscErrorCode SVDMatGetLocalSize(SVD,PetscInt*,PetscInt*);
1755 antodo 110
EXTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,Vec*,Vec,Vec*,PetscInt,PetscInt,PetscScalar*);
1521 slepc 111
 
1951 jroman 112
#endif