Actual source code: svdimpl.h
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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 PetscErrorCode SVDRegisterAll(void);
19: SLEPC_EXTERN PetscLogEvent SVD_SetUp,SVD_Solve;
21: typedef struct _SVDOps *SVDOps;
23: struct _SVDOps {
24: PetscErrorCode (*solve)(SVD);
25: PetscErrorCode (*setup)(SVD);
26: PetscErrorCode (*setfromoptions)(PetscOptionItems*,SVD);
27: PetscErrorCode (*publishoptions)(SVD);
28: PetscErrorCode (*destroy)(SVD);
29: PetscErrorCode (*reset)(SVD);
30: PetscErrorCode (*view)(SVD,PetscViewer);
31: };
33: /*
34: Maximum number of monitors you can run with a single SVD
35: */
36: #define MAXSVDMONITORS 5
38: typedef enum { SVD_STATE_INITIAL,
39: SVD_STATE_SETUP,
40: SVD_STATE_SOLVED,
41: SVD_STATE_VECTORS } SVDStateType;
43: /*
44: To check for unsupported features at SVDSetUp_XXX()
45: */
46: typedef enum { SVD_FEATURE_CONVERGENCE=16, /* convergence test selected by user */
47: SVD_FEATURE_STOPPING=32 /* stopping test */
48: } SVDFeatureType;
50: /*
51: Defines the SVD data structure.
52: */
53: struct _p_SVD {
54: PETSCHEADER(struct _SVDOps);
55: /*------------------------- User parameters ---------------------------*/
56: Mat OP; /* problem matrix */
57: PetscInt max_it; /* max iterations */
58: PetscInt nsv; /* number of requested values */
59: PetscInt ncv; /* basis size */
60: PetscInt mpd; /* maximum dimension of projected problem */
61: PetscInt nini,ninil; /* number of initial vecs (negative means not copied yet) */
62: PetscReal tol; /* tolerance */
63: SVDConv conv; /* convergence test */
64: SVDStop stop; /* stopping test */
65: SVDWhich which; /* which singular values are computed */
66: PetscBool impltrans; /* implicit transpose mode */
67: PetscBool trackall; /* whether all the residuals must be computed */
69: /*-------------- User-provided functions and contexts -----------------*/
70: PetscErrorCode (*converged)(SVD,PetscReal,PetscReal,PetscReal*,void*);
71: PetscErrorCode (*convergeduser)(SVD,PetscReal,PetscReal,PetscReal*,void*);
72: PetscErrorCode (*convergeddestroy)(void*);
73: PetscErrorCode (*stopping)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
74: PetscErrorCode (*stoppinguser)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
75: PetscErrorCode (*stoppingdestroy)(void*);
76: void *convergedctx;
77: void *stoppingctx;
78: PetscErrorCode (*monitor[MAXSVDMONITORS])(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
79: PetscErrorCode (*monitordestroy[MAXSVDMONITORS])(void**);
80: void *monitorcontext[MAXSVDMONITORS];
81: PetscInt numbermonitors;
83: /*----------------- Child objects and working data -------------------*/
84: DS ds; /* direct solver object */
85: BV U,V; /* left and right singular vectors */
86: SlepcSC sc; /* sorting criterion data */
87: Mat A; /* problem matrix (m>n) */
88: Mat AT; /* transposed matrix */
89: Vec *IS,*ISL; /* placeholder for references to user initial space */
90: PetscReal *sigma; /* singular values */
91: PetscInt *perm; /* permutation for singular value ordering */
92: PetscReal *errest; /* error estimates */
93: void *data; /* placeholder for solver-specific stuff */
95: /* ----------------------- Status variables -------------------------- */
96: SVDStateType state; /* initial -> setup -> solved -> vectors */
97: PetscInt nconv; /* number of converged values */
98: PetscInt its; /* iteration counter */
99: PetscBool leftbasis; /* if U is filled by the solver */
100: SVDConvergedReason reason;
101: };
103: /*
104: Macros to test valid SVD arguments
105: */
106: #if !defined(PETSC_USE_DEBUG)
108: #define SVDCheckSolved(h,arg) do {} while (0)
110: #else
112: #define SVDCheckSolved(h,arg) \
113: do { \
114: if ((h)->state<SVD_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call SVDSolve() first: Parameter #%d",arg); \
115: } while (0)
117: #endif
119: /* Check for unsupported features */
120: #define SVDCheckUnsupportedCondition(svd,mask,condition,msg) \
121: do { \
122: if (condition) { \
123: if (((mask) & SVD_FEATURE_CONVERGENCE) && (svd)->converged!=SVDConvergedRelative) SETERRQ2(PetscObjectComm((PetscObject)(svd)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(svd))->type_name,(msg)); \
124: if (((mask) & SVD_FEATURE_STOPPING) && (svd)->stopping!=SVDStoppingBasic) SETERRQ2(PetscObjectComm((PetscObject)(svd)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(svd))->type_name,(msg)); \
125: } \
126: } while (0)
127: #define SVDCheckUnsupported(svd,mask) SVDCheckUnsupportedCondition(svd,mask,PETSC_TRUE,"")
129: /* Check for ignored features */
130: #define SVDCheckIgnoredCondition(svd,mask,condition,msg) \
131: do { \
132: PetscErrorCode __ierr; \
133: if (condition) { \
134: if (((mask) & SVD_FEATURE_CONVERGENCE) && (svd)->converged!=SVDConvergedRelative) { __PetscInfo2((svd),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(svd))->type_name,(msg)); } \
135: if (((mask) & SVD_FEATURE_STOPPING) && (svd)->stopping!=SVDStoppingBasic) { __PetscInfo2((svd),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(svd))->type_name,(msg)); } \
136: } \
137: } while (0)
138: #define SVDCheckIgnored(svd,mask) SVDCheckIgnoredCondition(svd,mask,PETSC_TRUE,"")
140: PETSC_STATIC_INLINE PetscErrorCode SVDMatMult(SVD svd,PetscBool trans,Vec x,Vec y)
141: {
145: if (trans) {
146: if (svd->AT) {
147: MatMult(svd->AT,x,y);
148: } else {
149: #if defined(PETSC_USE_COMPLEX)
150: MatMultHermitianTranspose(svd->A,x,y);
151: #else
152: MatMultTranspose(svd->A,x,y);
153: #endif
154: }
155: } else {
156: if (svd->A) {
157: MatMult(svd->A,x,y);
158: } else {
159: #if defined(PETSC_USE_COMPLEX)
160: MatMultHermitianTranspose(svd->AT,x,y);
161: #else
162: MatMultTranspose(svd->AT,x,y);
163: #endif
164: }
165: }
166: return(0);
167: }
169: PETSC_STATIC_INLINE PetscErrorCode SVDMatCreateVecs(SVD svd,Vec *x,Vec *y)
170: {
174: if (svd->A) {
175: MatCreateVecs(svd->A,x,y);
176: } else {
177: MatCreateVecs(svd->AT,y,x);
178: }
179: return(0);
180: }
182: PETSC_STATIC_INLINE PetscErrorCode SVDMatCreateVecsEmpty(SVD svd,Vec *x,Vec *y)
183: {
187: if (svd->A) {
188: MatCreateVecsEmpty(svd->A,x,y);
189: } else {
190: MatCreateVecsEmpty(svd->AT,y,x);
191: }
192: return(0);
193: }
195: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetSize(SVD svd,PetscInt *m,PetscInt *n)
196: {
200: if (svd->A) {
201: MatGetSize(svd->A,m,n);
202: } else {
203: MatGetSize(svd->AT,n,m);
204: }
205: return(0);
206: }
208: PETSC_STATIC_INLINE PetscErrorCode SVDMatGetLocalSize(SVD svd,PetscInt *m,PetscInt *n)
209: {
213: if (svd->A) {
214: MatGetLocalSize(svd->A,m,n);
215: } else {
216: MatGetLocalSize(svd->AT,n,m);
217: }
218: return(0);
219: }
221: SLEPC_INTERN PetscErrorCode SVDTwoSideLanczos(SVD,PetscReal*,PetscReal*,BV,BV,PetscInt,PetscInt);
222: SLEPC_INTERN PetscErrorCode SVDSetDimensions_Default(SVD);
223: SLEPC_INTERN PetscErrorCode SVDComputeVectors(SVD);
225: #endif