Actual source code: slepcsvd.h
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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: */
10: /*
11: User interface for SLEPc's singular value solvers
12: */
16: #include <slepceps.h>
17: #include <slepcbv.h>
18: #include <slepcds.h>
20: PETSC_EXTERN PetscErrorCode SVDInitializePackage(void);
22: /*S
23: SVD - Abstract SLEPc object that manages all the singular value
24: problem solvers.
26: Level: beginner
28: .seealso: SVDCreate()
29: S*/
30: typedef struct _p_SVD* SVD;
32: /*J
33: SVDType - String with the name of a SLEPc singular value solver
35: Level: beginner
37: .seealso: SVDSetType(), SVD
38: J*/
39: typedef const char* SVDType;
40: #define SVDCROSS "cross"
41: #define SVDCYCLIC "cyclic"
42: #define SVDLAPACK "lapack"
43: #define SVDLANCZOS "lanczos"
44: #define SVDTRLANCZOS "trlanczos"
45: #define SVDPRIMME "primme"
47: /* Logging support */
48: PETSC_EXTERN PetscClassId SVD_CLASSID;
50: /*E
51: SVDWhich - Determines whether largest or smallest singular triplets
52: are to be computed
54: Level: intermediate
56: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
57: E*/
58: typedef enum { SVD_LARGEST,
59: SVD_SMALLEST } SVDWhich;
61: /*E
62: SVDErrorType - The error type used to assess accuracy of computed solutions
64: Level: intermediate
66: .seealso: SVDComputeError()
67: E*/
68: typedef enum { SVD_ERROR_ABSOLUTE,
69: SVD_ERROR_RELATIVE } SVDErrorType;
70: PETSC_EXTERN const char *SVDErrorTypes[];
72: /*E
73: SVDConv - Determines the convergence test
75: Level: intermediate
77: .seealso: SVDSetConvergenceTest(), SVDSetConvergenceTestFunction()
78: E*/
79: typedef enum { SVD_CONV_ABS,
80: SVD_CONV_REL,
81: SVD_CONV_USER } SVDConv;
83: /*E
84: SVDStop - Determines the stopping test
86: Level: advanced
88: .seealso: SVDSetStoppingTest(), SVDSetStoppingTestFunction()
89: E*/
90: typedef enum { SVD_STOP_BASIC,
91: SVD_STOP_USER } SVDStop;
93: /*E
94: SVDConvergedReason - Reason a singular value solver was said to
95: have converged or diverged
97: Level: intermediate
99: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
100: E*/
101: typedef enum {/* converged */
102: SVD_CONVERGED_TOL = 1,
103: SVD_CONVERGED_USER = 2,
104: /* diverged */
105: SVD_DIVERGED_ITS = -1,
106: SVD_DIVERGED_BREAKDOWN = -2,
107: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
108: PETSC_EXTERN const char *const*SVDConvergedReasons;
110: PETSC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
111: PETSC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
112: PETSC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
113: PETSC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
114: PETSC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
115: PETSC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
116: PETSC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
117: PETSC_EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
118: PETSC_EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
119: PETSC_EXTERN PetscErrorCode SVDSetInitialSpaces(SVD,PetscInt,Vec*,PetscInt,Vec*);
120: PETSC_DEPRECATED("Use SVDSetInitialSpaces()") PETSC_STATIC_INLINE PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt nr,Vec *isr) {return SVDSetInitialSpaces(svd,nr,isr,0,NULL);}
121: PETSC_DEPRECATED("Use SVDSetInitialSpaces()") PETSC_STATIC_INLINE PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt nl,Vec *isl) {return SVDSetInitialSpaces(svd,0,NULL,nl,isl);}
122: PETSC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
123: PETSC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
124: PETSC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
125: PETSC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
126: PETSC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
127: PETSC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
128: PETSC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
129: PETSC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
130: PETSC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
131: PETSC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
132: PETSC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
133: PETSC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
134: PETSC_EXTERN PetscErrorCode SVDSetUp(SVD);
135: PETSC_EXTERN PetscErrorCode SVDSolve(SVD);
136: PETSC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
137: PETSC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,PetscErrorCode (*)(SVD,PetscReal,PetscReal,PetscReal*,void*),void*,PetscErrorCode (*)(void*));
138: PETSC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
139: PETSC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
140: PETSC_EXTERN PetscErrorCode SVDConvergedAbsolute(SVD,PetscReal,PetscReal,PetscReal*,void*);
141: PETSC_EXTERN PetscErrorCode SVDConvergedRelative(SVD,PetscReal,PetscReal,PetscReal*,void*);
142: PETSC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
143: PETSC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
144: PETSC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
145: PETSC_EXTERN PetscErrorCode SVDStoppingBasic(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
146: PETSC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
147: PETSC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
148: PETSC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
149: PETSC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
150: PETSC_DEPRECATED("Use SVDComputeError()") PETSC_STATIC_INLINE PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *r) {return SVDComputeError(svd,i,SVD_ERROR_RELATIVE,r);}
151: PETSC_DEPRECATED("Use SVDComputeError() with SVD_ERROR_ABSOLUTE") PETSC_STATIC_INLINE PetscErrorCode SVDComputeResidualNorms(SVD svd,PetscInt i,PetscReal *r1,PETSC_UNUSED PetscReal *r2) {return SVDComputeError(svd,i,SVD_ERROR_ABSOLUTE,r1);}
152: PETSC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
153: PETSC_STATIC_INLINE PetscErrorCode SVDViewFromOptions(SVD svd,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)svd,obj,name);}
154: PETSC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
155: PETSC_DEPRECATED("Use SVDErrorView()") PETSC_STATIC_INLINE PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
156: PETSC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
157: PETSC_EXTERN PetscErrorCode SVDReasonView(SVD,PetscViewer);
158: PETSC_EXTERN PetscErrorCode SVDReasonViewFromOptions(SVD);
159: PETSC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
160: PETSC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
161: PETSC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
162: PETSC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
163: PETSC_EXTERN PetscErrorCode SVDDestroy(SVD*);
164: PETSC_EXTERN PetscErrorCode SVDReset(SVD);
166: PETSC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
167: PETSC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
168: PETSC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char*,const char*,const char*,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool);
169: PETSC_EXTERN PetscErrorCode SVDConvMonitorSetFromOptions(SVD,const char*,const char*,const char*,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,SlepcConvMonitor));
170: PETSC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
171: PETSC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
172: PETSC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
173: PETSC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
174: PETSC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,SlepcConvMonitor);
175: PETSC_EXTERN PetscErrorCode SVDMonitorLGCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
176: PETSC_EXTERN PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
177: PETSC_EXTERN PetscErrorCode SVDMonitorLGAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
179: PETSC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
180: PETSC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
182: PETSC_EXTERN PetscFunctionList SVDList;
183: PETSC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
185: PETSC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
187: /* --------- options specific to particular solvers -------- */
189: PETSC_EXTERN PetscErrorCode SVDCrossSetExplicitMatrix(SVD,PetscBool);
190: PETSC_EXTERN PetscErrorCode SVDCrossGetExplicitMatrix(SVD,PetscBool*);
191: PETSC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
192: PETSC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
194: PETSC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
195: PETSC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
196: PETSC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
197: PETSC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
199: PETSC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
200: PETSC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
202: PETSC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
203: PETSC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
205: /*E
206: SVDPRIMMEMethod - determines the SVD method selected in the PRIMME library
208: Level: advanced
210: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEGetMethod()
211: E*/
212: typedef enum { SVD_PRIMME_HYBRID=1,
213: SVD_PRIMME_NORMALEQUATIONS,
214: SVD_PRIMME_AUGMENTED } SVDPRIMMEMethod;
215: PETSC_EXTERN const char *SVDPRIMMEMethods[];
217: PETSC_EXTERN PetscErrorCode SVDPRIMMESetBlockSize(SVD,PetscInt);
218: PETSC_EXTERN PetscErrorCode SVDPRIMMEGetBlockSize(SVD,PetscInt*);
219: PETSC_EXTERN PetscErrorCode SVDPRIMMESetMethod(SVD,SVDPRIMMEMethod);
220: PETSC_EXTERN PetscErrorCode SVDPRIMMEGetMethod(SVD,SVDPRIMMEMethod*);
222: #endif