Actual source code: svdimpl.h

slepc-3.17.2 2022-08-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: #if !defined(SLEPCSVDIMPL_H)
 12: #define SLEPCSVDIMPL_H

 14: #include <slepcsvd.h>
 15: #include <slepc/private/slepcimpl.h>

 17: SLEPC_EXTERN PetscBool SVDRegisterAllCalled;
 18: SLEPC_EXTERN PetscBool SVDMonitorRegisterAllCalled;
 19: SLEPC_EXTERN PetscErrorCode SVDRegisterAll(void);
 20: SLEPC_EXTERN PetscErrorCode SVDMonitorRegisterAll(void);
 21: SLEPC_EXTERN PetscLogEvent SVD_SetUp,SVD_Solve;

 23: typedef struct _SVDOps *SVDOps;

 25: struct _SVDOps {
 26:   PetscErrorCode (*solve)(SVD);
 27:   PetscErrorCode (*solveg)(SVD);
 28:   PetscErrorCode (*setup)(SVD);
 29:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,SVD);
 30:   PetscErrorCode (*publishoptions)(SVD);
 31:   PetscErrorCode (*destroy)(SVD);
 32:   PetscErrorCode (*reset)(SVD);
 33:   PetscErrorCode (*view)(SVD,PetscViewer);
 34:   PetscErrorCode (*computevectors)(SVD);
 35: };

 37: /*
 38:      Maximum number of monitors you can run with a single SVD
 39: */
 40: #define MAXSVDMONITORS 5

 42: typedef enum { SVD_STATE_INITIAL,
 43:                SVD_STATE_SETUP,
 44:                SVD_STATE_SOLVED,
 45:                SVD_STATE_VECTORS } SVDStateType;

 47: /*
 48:    To check for unsupported features at SVDSetUp_XXX()
 49: */
 50: typedef enum { SVD_FEATURE_CONVERGENCE=16,  /* convergence test selected by user */
 51:                SVD_FEATURE_STOPPING=32      /* stopping test */
 52:              } SVDFeatureType;

 54: /*
 55:    Defines the SVD data structure.
 56: */
 57: struct _p_SVD {
 58:   PETSCHEADER(struct _SVDOps);
 59:   /*------------------------- User parameters ---------------------------*/
 60:   Mat            OP,OPb;           /* problem matrices */
 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:   SVDProblemType problem_type;     /* which kind of problem to be solved */
 71:   PetscBool      impltrans;        /* implicit transpose mode */
 72:   PetscBool      trackall;         /* whether all the residuals must be computed */

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

 88:   /*----------------- Child objects and working data -------------------*/
 89:   DS             ds;               /* direct solver object */
 90:   BV             U,V;              /* left and right singular vectors */
 91:   SlepcSC        sc;               /* sorting criterion data */
 92:   Mat            A,B;              /* problem matrices */
 93:   Mat            AT,BT;            /* transposed matrices */
 94:   Vec            *IS,*ISL;         /* placeholder for references to user initial space */
 95:   PetscReal      *sigma;           /* singular values */
 96:   PetscReal      *errest;          /* error estimates */
 97:   PetscInt       *perm;            /* permutation for singular value ordering */
 98:   PetscInt       nworkl,nworkr;    /* number of work vectors */
 99:   Vec            *workl,*workr;    /* work vectors */
100:   void           *data;            /* placeholder for solver-specific stuff */

102:   /* ----------------------- Status variables -------------------------- */
103:   SVDStateType   state;            /* initial -> setup -> solved -> vectors */
104:   PetscInt       nconv;            /* number of converged values */
105:   PetscInt       its;              /* iteration counter */
106:   PetscBool      leftbasis;        /* if U is filled by the solver */
107:   PetscBool      swapped;          /* the U and V bases have been swapped (M<N) */
108:   PetscBool      expltrans;        /* explicit transpose created */
109:   PetscReal      nrma,nrmb;        /* computed matrix norms */
110:   PetscBool      isgeneralized;
111:   SVDConvergedReason reason;
112: };

114: /*
115:     Macros to test valid SVD arguments
116: */
117: #if !defined(PETSC_USE_DEBUG)

119: #define SVDCheckSolved(h,arg) do {(void)(h);} while (0)

121: #else

123: #define SVDCheckSolved(h,arg) \
124:   do { \
126:   } while (0)

128: #endif

130: /*
131:     Macros to check settings at SVDSetUp()
132: */

134: /* SVDCheckStandard: the problem is not GSVD */
135: #define SVDCheckStandardCondition(svd,condition,msg) \
136:   do { \
137:     if (condition) { \
139:     } \
140:   } while (0)
141: #define SVDCheckStandard(svd) SVDCheckStandardCondition(svd,PETSC_TRUE,"")

143: /* Check for unsupported features */
144: #define SVDCheckUnsupportedCondition(svd,mask,condition,msg) \
145:   do { \
146:     if (condition) { \
149:     } \
150:   } while (0)
151: #define SVDCheckUnsupported(svd,mask) SVDCheckUnsupportedCondition(svd,mask,PETSC_TRUE,"")

153: /* Check for ignored features */
154: #define SVDCheckIgnoredCondition(svd,mask,condition,msg) \
155:   do { \
156:     if (condition) { \
157:       if (((mask) & SVD_FEATURE_CONVERGENCE) && (svd)->converged!=SVDConvergedRelative) PetscInfo((svd),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(svd))->type_name,(msg)); \
158:       if (((mask) & SVD_FEATURE_STOPPING) && (svd)->stopping!=SVDStoppingBasic) PetscInfo((svd),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(svd))->type_name,(msg)); \
159:     } \
160:   } while (0)
161: #define SVDCheckIgnored(svd,mask) SVDCheckIgnoredCondition(svd,mask,PETSC_TRUE,"")

163: /*
164:   SVD_KSPSetOperators - Sets the KSP matrices
165: */
166: static inline PetscErrorCode SVD_KSPSetOperators(KSP ksp,Mat A,Mat B)
167: {
168:   const char     *prefix;

170:   KSPSetOperators(ksp,A,B);
171:   MatGetOptionsPrefix(B,&prefix);
172:   if (!prefix) {
173:     /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS)
174:        only applies if the Mat has no user-defined prefix */
175:     KSPGetOptionsPrefix(ksp,&prefix);
176:     MatSetOptionsPrefix(B,prefix);
177:   }
178:   PetscFunctionReturn(0);
179: }

181: /*
182:    Create the template vector for the left basis in GSVD, as in
183:    MatCreateVecsEmpty(Z,NULL,&t) for Z=[A;B] without forming Z.
184:  */
185: static inline PetscErrorCode SVDCreateLeftTemplate(SVD svd,Vec *t)
186: {
187:   PetscInt       M,P,m,p;
188:   Vec            v1,v2;
189:   VecType        vec_type;

191:   MatCreateVecsEmpty(svd->OP,NULL,&v1);
192:   VecGetSize(v1,&M);
193:   VecGetLocalSize(v1,&m);
194:   VecGetType(v1,&vec_type);
195:   MatCreateVecsEmpty(svd->OPb,NULL,&v2);
196:   VecGetSize(v2,&P);
197:   VecGetLocalSize(v2,&p);
198:   VecCreate(PetscObjectComm((PetscObject)(v1)),t);
199:   VecSetType(*t,vec_type);
200:   VecSetSizes(*t,m+p,M+P);
201:   VecSetUp(*t);
202:   VecDestroy(&v1);
203:   VecDestroy(&v2);
204:   PetscFunctionReturn(0);
205: }

207: SLEPC_INTERN PetscErrorCode SVDKrylovConvergence(SVD,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
208: SLEPC_INTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,BV,BV,PetscInt,PetscInt*,PetscBool*);
209: SLEPC_INTERN PetscErrorCode SVDSetDimensions_Default(SVD);
210: SLEPC_INTERN PetscErrorCode SVDComputeVectors(SVD);
211: SLEPC_INTERN PetscErrorCode SVDComputeVectors_Left(SVD);

213: #endif