Subversion Repositories slepc-dev

Rev

Go to most recent revision | 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
2116 eromero 4
   Copyright (c) 2002-2010, 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
 
2283 jroman 25
#include <slepcsvd.h>
26
#include <slepcip.h>
1521 slepc 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;
2216 jroman 77
  PetscBool        trackall;
1521 slepc 78
 
79
  PetscErrorCode  (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
80
  PetscErrorCode  (*monitordestroy[MAXSVDMONITORS])(void*);
81
  void            *monitorcontext[MAXSVDMONITORS];
82
  PetscInt        numbermonitors;
83
 
84
  PetscInt        matvecs;
85
};
86
 
2240 jroman 87
extern PetscErrorCode SVDRegisterAll(const char *);
88
extern PetscErrorCode SVDInitializePackage(const char*);
89
extern PetscErrorCode SVDFinalizePackage(void);
1521 slepc 90
 
91
#define SVDMonitor(svd,it,nconv,sigma,errest,nest) \
92
        { PetscErrorCode _ierr; PetscInt _i,_im = svd->numbermonitors; \
93
          for ( _i=0; _i<_im; _i++ ) {\
94
            _ierr=(*svd->monitor[_i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[_i]);\
95
            CHKERRQ(_ierr); \
96
          } \
97
        }
98
 
1951 jroman 99
/* context for SVDMonitorConverged */
100
typedef struct {
101
  PetscViewerASCIIMonitor viewer;
102
  PetscInt oldnconv;
103
} SVDMONITOR_CONV;
2240 jroman 104
extern PetscErrorCode SVDMonitorDestroy_Converged(SVDMONITOR_CONV*);
1951 jroman 105
 
2240 jroman 106
extern PetscErrorCode SVDDestroy_Default(SVD);
107
extern PetscErrorCode SVDMatMult(SVD,PetscBool,Vec,Vec);
108
extern PetscErrorCode SVDMatGetVecs(SVD,Vec*,Vec*);
109
extern PetscErrorCode SVDMatGetSize(SVD,PetscInt*,PetscInt*);
110
extern PetscErrorCode SVDMatGetLocalSize(SVD,PetscInt*,PetscInt*);
111
extern PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,Vec*,Vec,Vec*,PetscInt,PetscInt,PetscScalar*);
1521 slepc 112
 
1951 jroman 113
#endif