Actual source code: slepcsvd.h

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    User interface for SLEPc's singular value solvers.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 26: #include <slepceps.h>
 27: #include <slepcbv.h>
 28: #include <slepcds.h>

 30: PETSC_EXTERN PetscErrorCode SVDInitializePackage(void);

 32: /*S
 33:      SVD - Abstract SLEPc object that manages all the singular value
 34:      problem solvers.

 36:    Level: beginner

 38: .seealso:  SVDCreate()
 39: S*/
 40: typedef struct _p_SVD* SVD;

 42: /*J
 43:     SVDType - String with the name of a SLEPc singular value solver

 45:    Level: beginner

 47: .seealso: SVDSetType(), SVD
 48: J*/
 49: typedef const char* SVDType;
 50: #define SVDCROSS       "cross"
 51: #define SVDCYCLIC      "cyclic"
 52: #define SVDLAPACK      "lapack"
 53: #define SVDLANCZOS     "lanczos"
 54: #define SVDTRLANCZOS   "trlanczos"

 56: /* Logging support */
 57: PETSC_EXTERN PetscClassId SVD_CLASSID;

 59: /*E
 60:     SVDWhich - Determines whether largest or smallest singular triplets
 61:     are to be computed

 63:     Level: intermediate

 65: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
 66: E*/
 67: typedef enum { SVD_LARGEST,
 68:                SVD_SMALLEST } SVDWhich;

 70: /*E
 71:     SVDErrorType - The error type used to assess accuracy of computed solutions

 73:     Level: intermediate

 75: .seealso: SVDComputeError()
 76: E*/
 77: typedef enum { SVD_ERROR_ABSOLUTE,
 78:                SVD_ERROR_RELATIVE } SVDErrorType;
 79: PETSC_EXTERN const char *SVDErrorTypes[];

 81: /*E
 82:     SVDConv - Determines the convergence test

 84:     Level: intermediate

 86: .seealso: SVDSetConvergenceTest(), SVDSetConvergenceTestFunction()
 87: E*/
 88: typedef enum { SVD_CONV_ABS,
 89:                SVD_CONV_REL,
 90:                SVD_CONV_USER } SVDConv;

 92: /*E
 93:     SVDStop - Determines the stopping test

 95:     Level: advanced

 97: .seealso: SVDSetStoppingTest(), SVDSetStoppingTestFunction()
 98: E*/
 99: typedef enum { SVD_STOP_BASIC,
100:                SVD_STOP_USER } SVDStop;

102: /*E
103:     SVDConvergedReason - Reason a singular value solver was said to
104:          have converged or diverged

106:    Level: intermediate

108: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
109: E*/
110: typedef enum {/* converged */
111:               SVD_CONVERGED_TOL                =  1,
112:               SVD_CONVERGED_USER               =  2,
113:               /* diverged */
114:               SVD_DIVERGED_ITS                 = -1,
115:               SVD_DIVERGED_BREAKDOWN           = -2,
116:               SVD_CONVERGED_ITERATING          =  0 } SVDConvergedReason;
117: PETSC_EXTERN const char *const*SVDConvergedReasons;

119: PETSC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
120: PETSC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
121: PETSC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
122: PETSC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
123: PETSC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
124: PETSC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
125: PETSC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
126: PETSC_EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
127: PETSC_EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
128: PETSC_EXTERN PetscErrorCode SVDSetInitialSpace(SVD,PetscInt,Vec*);
129: PETSC_EXTERN PetscErrorCode SVDSetInitialSpaceLeft(SVD,PetscInt,Vec*);
130: PETSC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
131: PETSC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
132: PETSC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
133: PETSC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
134: PETSC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
135: PETSC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
136: PETSC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
137: PETSC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
138: PETSC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
139: PETSC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
140: PETSC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
141: PETSC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
142: PETSC_EXTERN PetscErrorCode SVDSetUp(SVD);
143: PETSC_EXTERN PetscErrorCode SVDSolve(SVD);
144: PETSC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
145: PETSC_EXTERN PetscErrorCode SVDSetConvergenceTestFunction(SVD,PetscErrorCode (*)(SVD,PetscReal,PetscReal,PetscReal*,void*),void*,PetscErrorCode (*)(void*));
146: PETSC_EXTERN PetscErrorCode SVDSetConvergenceTest(SVD,SVDConv);
147: PETSC_EXTERN PetscErrorCode SVDGetConvergenceTest(SVD,SVDConv*);
148: PETSC_EXTERN PetscErrorCode SVDConvergedAbsolute(SVD,PetscReal,PetscReal,PetscReal*,void*);
149: PETSC_EXTERN PetscErrorCode SVDConvergedRelative(SVD,PetscReal,PetscReal,PetscReal*,void*);
150: PETSC_EXTERN PetscErrorCode SVDSetStoppingTestFunction(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
151: PETSC_EXTERN PetscErrorCode SVDSetStoppingTest(SVD,SVDStop);
152: PETSC_EXTERN PetscErrorCode SVDGetStoppingTest(SVD,SVDStop*);
153: PETSC_EXTERN PetscErrorCode SVDStoppingBasic(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*);
154: PETSC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
155: PETSC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
156: PETSC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
157: PETSC_EXTERN PetscErrorCode SVDComputeError(SVD,PetscInt,SVDErrorType,PetscReal*);
158: PETSC_DEPRECATED("Use SVDComputeError()") PETSC_STATIC_INLINE PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *r) {return SVDComputeError(svd,i,SVD_ERROR_RELATIVE,r);}
159: PETSC_DEPRECATED("Use SVDComputeError() with SVD_ERROR_ABSOLUTE") PETSC_STATIC_INLINE PetscErrorCode SVDComputeResidualNorms(SVD svd,PetscInt i,PetscReal *r1,PetscReal *r2) {return SVDComputeError(svd,i,SVD_ERROR_ABSOLUTE,r1);}
160: PETSC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
161: PETSC_STATIC_INLINE PetscErrorCode SVDViewFromOptions(SVD svd,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)svd,obj,name);}
162: PETSC_EXTERN PetscErrorCode SVDErrorView(SVD,SVDErrorType,PetscViewer);
163: PETSC_DEPRECATED("Use SVDErrorView()") PETSC_STATIC_INLINE PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer v) {return SVDErrorView(svd,SVD_ERROR_RELATIVE,v);}
164: PETSC_EXTERN PetscErrorCode SVDErrorViewFromOptions(SVD);
165: PETSC_EXTERN PetscErrorCode SVDReasonView(SVD,PetscViewer);
166: PETSC_EXTERN PetscErrorCode SVDReasonViewFromOptions(SVD);
167: PETSC_EXTERN PetscErrorCode SVDValuesView(SVD,PetscViewer);
168: PETSC_EXTERN PetscErrorCode SVDValuesViewFromOptions(SVD);
169: PETSC_EXTERN PetscErrorCode SVDVectorsView(SVD,PetscViewer);
170: PETSC_EXTERN PetscErrorCode SVDVectorsViewFromOptions(SVD);
171: PETSC_EXTERN PetscErrorCode SVDDestroy(SVD*);
172: PETSC_EXTERN PetscErrorCode SVDReset(SVD);

174: PETSC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
175: PETSC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
176: PETSC_EXTERN PetscErrorCode SVDMonitorSetFromOptions(SVD,const char*,const char*,const char*,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool);
177: PETSC_EXTERN PetscErrorCode SVDConvMonitorSetFromOptions(SVD,const char*,const char*,const char*,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,SlepcConvMonitor));
178: PETSC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
179: PETSC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
180: PETSC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
181: PETSC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*);
182: PETSC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,SlepcConvMonitor);
183: PETSC_EXTERN PetscErrorCode SVDMonitorLGCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
184: PETSC_EXTERN PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
185: PETSC_EXTERN PetscErrorCode SVDMonitorLGAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);

187: PETSC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
188: PETSC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);

190: PETSC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
191: PETSC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);

193: PETSC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
194: PETSC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
195: PETSC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
196: PETSC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);

198: PETSC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
199: PETSC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);

201: PETSC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
202: PETSC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);

204: PETSC_EXTERN PetscFunctionList SVDList;
205: PETSC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));

207: PETSC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);

209: #endif