Actual source code: svdimpl.h

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, 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:   /*------------------------- User parameters ---------------------------*/
 53:   Mat              OP;          /* problem matrix */
 54:   PetscInt         max_it;      /* max iterations */
 55:   PetscInt         nsv;         /* number of requested values */
 56:   PetscInt         ncv;         /* basis size */
 57:   PetscInt         mpd;         /* maximum dimension of projected problem */
 58:   PetscInt         nini,ninil;  /* number of initial vectors (negative means not copied yet) */
 59:   PetscReal        tol;         /* tolerance */
 60:   SVDWhich         which;       /* which singular values are computed */
 61:   PetscBool        impltrans;   /* implicit transpose mode */
 62:   PetscBool        trackall;    /* whether all the residuals must be computed */

 64:   /*-------------- User-provided functions and contexts -----------------*/
 65:   PetscErrorCode   (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
 66:   PetscErrorCode   (*monitordestroy[MAXSVDMONITORS])(void**);
 67:   void             *monitorcontext[MAXSVDMONITORS];
 68:   PetscInt         numbermonitors;

 70:   /*----------------- Child objects and working data -------------------*/
 71:   DS               ds;          /* direct solver object */
 72:   BV               U,V;         /* left and right singular vectors */
 73:   PetscRandom      rand;        /* random number generator */
 74:   SlepcSC          sc;          /* sorting criterion data */
 75:   Mat              A;           /* problem matrix (m>n) */
 76:   Mat              AT;          /* transposed matrix */
 77:   Vec              *IS,*ISL;    /* placeholder for references to user-provided initial space */
 78:   PetscReal        *sigma;      /* singular values */
 79:   PetscInt         *perm;       /* permutation for singular value ordering */
 80:   PetscReal        *errest;     /* error estimates */
 81:   void             *data;       /* placeholder for solver-specific stuff */

 83:   /* ----------------------- Status variables -------------------------- */
 84:   PetscInt         nconv;       /* number of converged values */
 85:   PetscInt         its;         /* iteration counter */
 86:   PetscBool        leftbasis;   /* if U is filled by the solver */
 87:   PetscBool        lvecsavail;  /* if U contains left singular vectors */
 88:   PetscInt         setupcalled;
 89:   SVDConvergedReason reason;
 90: };

 94: PETSC_STATIC_INLINE PetscErrorCode SVDMatMult(SVD svd,PetscBool trans,Vec x,Vec y)
 95: {

 99:   if (trans) {
100:     if (svd->AT) {
101:       MatMult(svd->AT,x,y);
102:     } else {
103: #if defined(PETSC_USE_COMPLEX)
104:       MatMultHermitianTranspose(svd->A,x,y);
105: #else
106:       MatMultTranspose(svd->A,x,y);
107: #endif
108:     }
109:   } else {
110:     if (svd->A) {
111:       MatMult(svd->A,x,y);
112:     } else {
113: #if defined(PETSC_USE_COMPLEX)
114:       MatMultHermitianTranspose(svd->AT,x,y);
115: #else
116:       MatMultTranspose(svd->AT,x,y);
117: #endif
118:     }
119:   }
120:   return(0);
121: }

125: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetVecs(SVD svd,Vec *x,Vec *y)
126: {

130:   if (svd->A) {
131:     MatGetVecs(svd->A,x,y);
132:   } else {
133:     MatGetVecs(svd->AT,y,x);
134:   }
135:   return(0);
136: }

140: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetSize(SVD svd,PetscInt *m,PetscInt *n)
141: {

145:   if (svd->A) {
146:     MatGetSize(svd->A,m,n);
147:   } else {
148:     MatGetSize(svd->AT,n,m);
149:   }
150:   return(0);
151: }

155: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetLocalSize(SVD svd,PetscInt *m,PetscInt *n)
156: {

160:   if (svd->A) {
161:     MatGetLocalSize(svd->A,m,n);
162:   } else {
163:     MatGetLocalSize(svd->AT,n,m);
164:   }
165:   return(0);
166: } 

168: PETSC_INTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,BV,BV,PetscInt,PetscInt);
169: PETSC_INTERN PetscErrorCode SVDSetDimensions_Default(SVD);

171: #endif