Actual source code: svdimpl.h

slepc-3.19.2 2023-09-05
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: /* SUBMANSEC = SVD */

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

 25: typedef struct _SVDOps *SVDOps;

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

 41: /*
 42:      Maximum number of monitors you can run with a single SVD
 43: */
 44: #define MAXSVDMONITORS 5

 46: typedef enum { SVD_STATE_INITIAL,
 47:                SVD_STATE_SETUP,
 48:                SVD_STATE_SOLVED,
 49:                SVD_STATE_VECTORS } SVDStateType;

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

 58: /*
 59:    Defines the SVD data structure.
 60: */
 61: struct _p_SVD {
 62:   PETSCHEADER(struct _SVDOps);
 63:   /*------------------------- User parameters ---------------------------*/
 64:   Mat            OP,OPb;           /* problem matrices */
 65:   Vec            omega;            /* signature for hyperbolic problems */
 66:   PetscInt       max_it;           /* max iterations */
 67:   PetscInt       nsv;              /* number of requested values */
 68:   PetscInt       ncv;              /* basis size */
 69:   PetscInt       mpd;              /* maximum dimension of projected problem */
 70:   PetscInt       nini,ninil;       /* number of initial vecs (negative means not copied yet) */
 71:   PetscReal      tol;              /* tolerance */
 72:   SVDConv        conv;             /* convergence test */
 73:   SVDStop        stop;             /* stopping test */
 74:   SVDWhich       which;            /* which singular values are computed */
 75:   SVDProblemType problem_type;     /* which kind of problem to be solved */
 76:   PetscBool      impltrans;        /* implicit transpose mode */
 77:   PetscBool      trackall;         /* whether all the residuals must be computed */

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

 93:   /*----------------- Child objects and working data -------------------*/
 94:   DS             ds;               /* direct solver object */
 95:   BV             U,V;              /* left and right singular vectors */
 96:   SlepcSC        sc;               /* sorting criterion data */
 97:   Mat            A,B;              /* problem matrices */
 98:   Mat            AT,BT;            /* transposed matrices */
 99:   Vec            *IS,*ISL;         /* placeholder for references to user initial space */
100:   PetscReal      *sigma;           /* singular values */
101:   PetscReal      *errest;          /* error estimates */
102:   PetscReal      *sign;            /* +-1 for each singular value in hyperbolic problems=U'*Omega*U */
103:   PetscInt       *perm;            /* permutation for singular value ordering */
104:   PetscInt       nworkl,nworkr;    /* number of work vectors */
105:   Vec            *workl,*workr;    /* work vectors */
106:   void           *data;            /* placeholder for solver-specific stuff */

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

121: /*
122:     Macros to test valid SVD arguments
123: */
124: #if !defined(PETSC_USE_DEBUG)

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

128: #else

130: #define SVDCheckSolved(h,arg) \
131:   do { \
132:     PetscCheck((h)->state>=SVD_STATE_SOLVED,PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call SVDSolve() first: Parameter #%d",arg); \
133:   } while (0)

135: #endif

137: /*
138:     Macros to check settings at SVDSetUp()
139: */

141: /* SVDCheckStandard: the problem is not GSVD */
142: #define SVDCheckStandardCondition(svd,condition,msg) \
143:   do { \
144:     if (condition) { \
145:       PetscCheck(!(svd)->isgeneralized,PetscObjectComm((PetscObject)(svd)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for generalized problems",((PetscObject)(svd))->type_name,(msg)); \
146:     } \
147:   } while (0)
148: #define SVDCheckStandard(svd) SVDCheckStandardCondition(svd,PETSC_TRUE,"")

150: /* SVDCheckDefinite: the problem is not hyperbolic */
151: #define SVDCheckDefiniteCondition(svd,condition,msg) \
152:   do { \
153:     if (condition) { \
154:       PetscCheck(!(svd)->ishyperbolic,PetscObjectComm((PetscObject)(svd)),PETSC_ERR_SUP,"The solver '%s'%s cannot be used for hyperbolic problems",((PetscObject)(svd))->type_name,(msg)); \
155:     } \
156:   } while (0)
157: #define SVDCheckDefinite(svd) SVDCheckDefiniteCondition(svd,PETSC_TRUE,"")

159: /* Check for unsupported features */
160: #define SVDCheckUnsupportedCondition(svd,mask,condition,msg) \
161:   do { \
162:     if (condition) { \
163:       PetscCheck(!((mask) & SVD_FEATURE_CONVERGENCE) || (svd)->converged==SVDConvergedRelative,PetscObjectComm((PetscObject)(svd)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(svd))->type_name,(msg)); \
164:       PetscCheck(!((mask) & SVD_FEATURE_STOPPING) || (svd)->stopping==SVDStoppingBasic,PetscObjectComm((PetscObject)(svd)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(svd))->type_name,(msg)); \
165:     } \
166:   } while (0)
167: #define SVDCheckUnsupported(svd,mask) SVDCheckUnsupportedCondition(svd,mask,PETSC_TRUE,"")

169: /* Check for ignored features */
170: #define SVDCheckIgnoredCondition(svd,mask,condition,msg) \
171:   do { \
172:     if (condition) { \
173:       if (((mask) & SVD_FEATURE_CONVERGENCE) && (svd)->converged!=SVDConvergedRelative) PetscCall(PetscInfo((svd),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(svd))->type_name,(msg))); \
174:       if (((mask) & SVD_FEATURE_STOPPING) && (svd)->stopping!=SVDStoppingBasic) PetscCall(PetscInfo((svd),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(svd))->type_name,(msg))); \
175:     } \
176:   } while (0)
177: #define SVDCheckIgnored(svd,mask) SVDCheckIgnoredCondition(svd,mask,PETSC_TRUE,"")

179: /*
180:   SVD_KSPSetOperators - Sets the KSP matrices
181: */
182: static inline PetscErrorCode SVD_KSPSetOperators(KSP ksp,Mat A,Mat B)
183: {
184:   const char     *prefix;

186:   PetscFunctionBegin;
187:   PetscCall(KSPSetOperators(ksp,A,B));
188:   PetscCall(MatGetOptionsPrefix(B,&prefix));
189:   if (!prefix) {
190:     /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS)
191:        only applies if the Mat has no user-defined prefix */
192:     PetscCall(KSPGetOptionsPrefix(ksp,&prefix));
193:     PetscCall(MatSetOptionsPrefix(B,prefix));
194:   }
195:   PetscFunctionReturn(PETSC_SUCCESS);
196: }

198: /*
199:    Create the template vector for the left basis in GSVD, as in
200:    MatCreateVecsEmpty(Z,NULL,&t) for Z=[A;B] without forming Z.
201:  */
202: static inline PetscErrorCode SVDCreateLeftTemplate(SVD svd,Vec *t)
203: {
204:   PetscInt       M,P,m,p;
205:   Vec            v1,v2;
206:   VecType        vec_type;

208:   PetscFunctionBegin;
209:   PetscCall(MatCreateVecsEmpty(svd->OP,NULL,&v1));
210:   PetscCall(VecGetSize(v1,&M));
211:   PetscCall(VecGetLocalSize(v1,&m));
212:   PetscCall(VecGetType(v1,&vec_type));
213:   PetscCall(MatCreateVecsEmpty(svd->OPb,NULL,&v2));
214:   PetscCall(VecGetSize(v2,&P));
215:   PetscCall(VecGetLocalSize(v2,&p));
216:   PetscCall(VecCreate(PetscObjectComm((PetscObject)(v1)),t));
217:   PetscCall(VecSetType(*t,vec_type));
218:   PetscCall(VecSetSizes(*t,m+p,M+P));
219:   PetscCall(VecSetUp(*t));
220:   PetscCall(VecDestroy(&v1));
221:   PetscCall(VecDestroy(&v2));
222:   PetscFunctionReturn(PETSC_SUCCESS);
223: }

225: SLEPC_INTERN PetscErrorCode SVDKrylovConvergence(SVD,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
226: SLEPC_INTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,BV,BV,PetscInt,PetscInt*,PetscBool*);
227: SLEPC_INTERN PetscErrorCode SVDSetDimensions_Default(SVD);
228: SLEPC_INTERN PetscErrorCode SVDComputeVectors(SVD);
229: SLEPC_INTERN PetscErrorCode SVDComputeVectors_Left(SVD);

231: #endif