Actual source code: svdimpl.h

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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 PetscBool SVDRegisterAllCalled;
 29: PETSC_EXTERN PetscErrorCode SVDRegisterAll(void);
 30: PETSC_EXTERN PetscLogEvent SVD_SetUp,SVD_Solve;

 32: typedef struct _SVDOps *SVDOps;

 34: struct _SVDOps {
 35:   PetscErrorCode (*solve)(SVD);
 36:   PetscErrorCode (*setup)(SVD);
 37:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,SVD);
 38:   PetscErrorCode (*publishoptions)(SVD);
 39:   PetscErrorCode (*destroy)(SVD);
 40:   PetscErrorCode (*reset)(SVD);
 41:   PetscErrorCode (*view)(SVD,PetscViewer);
 42: };

 44: /*
 45:      Maximum number of monitors you can run with a single SVD
 46: */
 47: #define MAXSVDMONITORS 5

 49: typedef enum { SVD_STATE_INITIAL,
 50:                SVD_STATE_SETUP,
 51:                SVD_STATE_SOLVED,
 52:                SVD_STATE_VECTORS } SVDStateType;

 54: /*
 55:    Defines the SVD data structure.
 56: */
 57: struct _p_SVD {
 58:   PETSCHEADER(struct _SVDOps);
 59:   /*------------------------- User parameters ---------------------------*/
 60:   Mat            OP;               /* problem matrix */
 61:   PetscInt       max_it;           /* max iterations */
 62:   PetscInt       nsv;              /* number of requested values */
 63:   PetscInt       ncv;              /* basis size */
 64:   PetscInt       mpd;              /* maximum dimension of projected problem */
 65:   PetscInt       nini,ninil;       /* number of initial vecs (negative means not copied yet) */
 66:   PetscReal      tol;              /* tolerance */
 67:   SVDConv        conv;             /* convergence test */
 68:   SVDStop        stop;             /* stopping test */
 69:   SVDWhich       which;            /* which singular values are computed */
 70:   PetscBool      impltrans;        /* implicit transpose mode */
 71:   PetscBool      trackall;         /* whether all the residuals must be computed */

 73:   /*-------------- User-provided functions and contexts -----------------*/
 74:   PetscErrorCode (*converged)(SVD,PetscReal,PetscReal,PetscReal*,void*);
 75:   PetscErrorCode (*convergeddestroy)(void*);
 76:   PetscErrorCode (*stopping)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
 77:   PetscErrorCode (*stoppingdestroy)(void*);
 78:   void           *convergedctx;
 79:   void           *stoppingctx;
 80:   PetscErrorCode (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
 81:   PetscErrorCode (*monitordestroy[MAXSVDMONITORS])(void**);
 82:   void           *monitorcontext[MAXSVDMONITORS];
 83:   PetscInt       numbermonitors;

 85:   /*----------------- Child objects and working data -------------------*/
 86:   DS             ds;               /* direct solver object */
 87:   BV             U,V;              /* left and right singular vectors */
 88:   SlepcSC        sc;               /* sorting criterion data */
 89:   Mat            A;                /* problem matrix (m>n) */
 90:   Mat            AT;               /* transposed matrix */
 91:   Vec            *IS,*ISL;         /* placeholder for references to user initial space */
 92:   PetscReal      *sigma;           /* singular values */
 93:   PetscInt       *perm;            /* permutation for singular value ordering */
 94:   PetscReal      *errest;          /* error estimates */
 95:   void           *data;            /* placeholder for solver-specific stuff */

 97:   /* ----------------------- Status variables -------------------------- */
 98:   SVDStateType   state;            /* initial -> setup -> solved -> vectors */
 99:   PetscInt       nconv;            /* number of converged values */
100:   PetscInt       its;              /* iteration counter */
101:   PetscBool      leftbasis;        /* if U is filled by the solver */
102:   SVDConvergedReason reason;
103: };

105: /*
106:     Macros to test valid SVD arguments
107: */
108: #if !defined(PETSC_USE_DEBUG)

110: #define SVDCheckSolved(h,arg) do {} while (0)

112: #else

114: #define SVDCheckSolved(h,arg) \
115:   do { \
116:     if (h->state<SVD_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call SVDSolve() first: Parameter #%d",arg); \
117:   } while (0)

119: #endif

123: PETSC_STATIC_INLINE PetscErrorCode SVDMatMult(SVD svd,PetscBool trans,Vec x,Vec y)
124: {

128:   if (trans) {
129:     if (svd->AT) {
130:       MatMult(svd->AT,x,y);
131:     } else {
132: #if defined(PETSC_USE_COMPLEX)
133:       MatMultHermitianTranspose(svd->A,x,y);
134: #else
135:       MatMultTranspose(svd->A,x,y);
136: #endif
137:     }
138:   } else {
139:     if (svd->A) {
140:       MatMult(svd->A,x,y);
141:     } else {
142: #if defined(PETSC_USE_COMPLEX)
143:       MatMultHermitianTranspose(svd->AT,x,y);
144: #else
145:       MatMultTranspose(svd->AT,x,y);
146: #endif
147:     }
148:   }
149:   return(0);
150: }

154: PETSC_STATIC_INLINE PetscErrorCode SVDMatCreateVecs(SVD svd,Vec *x,Vec *y)
155: {

159:   if (svd->A) {
160:     MatCreateVecs(svd->A,x,y);
161:   } else {
162:     MatCreateVecs(svd->AT,y,x);
163:   }
164:   return(0);
165: }

169: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetSize(SVD svd,PetscInt *m,PetscInt *n)
170: {

174:   if (svd->A) {
175:     MatGetSize(svd->A,m,n);
176:   } else {
177:     MatGetSize(svd->AT,n,m);
178:   }
179:   return(0);
180: }

184: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetLocalSize(SVD svd,PetscInt *m,PetscInt *n)
185: {

189:   if (svd->A) {
190:     MatGetLocalSize(svd->A,m,n);
191:   } else {
192:     MatGetLocalSize(svd->AT,n,m);
193:   }
194:   return(0);
195: } 

197: PETSC_INTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,BV,BV,PetscInt,PetscInt);
198: PETSC_INTERN PetscErrorCode SVDSetDimensions_Default(SVD);
199: PETSC_INTERN PetscErrorCode SVDComputeVectors(SVD);

201: #endif