Actual source code: slepcds.h

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

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

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

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

 24: #include <slepcsc.h>
 25: #include <slepcfn.h>

 27: #define DS_MAX_SOLVE 6

 29: PETSC_EXTERN PetscErrorCode DSInitializePackage(void);
 30: /*S
 31:     DS - Direct solver (or dense system), to represent low-dimensional
 32:     eigenproblems that must be solved within iterative solvers. This is an
 33:     auxiliary object and is not normally needed by application programmers.

 35:     Level: beginner

 37: .seealso:  DSCreate()
 38: S*/
 39: typedef struct _p_DS* DS;

 41: /*J
 42:     DSType - String with the name of the type of direct solver. Roughly,
 43:     there are as many types as problem types are available within SLEPc.

 45:     Level: advanced

 47: .seealso: DSSetType(), DS
 48: J*/
 49: typedef const char* DSType;
 50: #define DSHEP             "hep"
 51: #define DSNHEP            "nhep"
 52: #define DSGHEP            "ghep"
 53: #define DSGHIEP           "ghiep"
 54: #define DSGNHEP           "gnhep"
 55: #define DSSVD             "svd"
 56: #define DSPEP             "pep"
 57: #define DSNEP             "nep"

 59: /* Logging support */
 60: PETSC_EXTERN PetscClassId DS_CLASSID;

 62: /*E
 63:     DSStateType - Indicates in which state the direct solver is

 65:     Level: advanced

 67: .seealso: DSSetState()
 68: E*/
 69: typedef enum { DS_STATE_RAW,
 70:                DS_STATE_INTERMEDIATE,
 71:                DS_STATE_CONDENSED,
 72:                DS_STATE_TRUNCATED } DSStateType;

 74: /*E
 75:     DSMatType - Used to refer to one of the matrices stored internally in DS

 77:     Notes:
 78:     The matrices preferently refer to
 79: +   DS_MAT_A  - first matrix of eigenproblem/singular value problem
 80: .   DS_MAT_B  - second matrix of a generalized eigenproblem
 81: .   DS_MAT_C  - third matrix of a quadratic eigenproblem
 82: .   DS_MAT_T  - tridiagonal matrix
 83: .   DS_MAT_D  - diagonal matrix
 84: .   DS_MAT_Q  - orthogonal matrix of (right) Schur vectors
 85: .   DS_MAT_Z  - orthogonal matrix of left Schur vectors
 86: .   DS_MAT_X  - right eigenvectors
 87: .   DS_MAT_Y  - left eigenvectors
 88: .   DS_MAT_U  - left singular vectors
 89: .   DS_MAT_VT - right singular vectors
 90: .   DS_MAT_W  - workspace matrix
 91: -   DS_MAT_Ex - extra matrices (x=0,..,9)

 93:     All matrices can have space to hold ld x ld elements, except for
 94:     DS_MAT_T that has space for 3 x ld elements (ld = leading dimension)
 95:     and DS_MAT_D that has space for just ld elements.

 97:     In DSPEP problems, matrices A, B, W can have space for d*ld x d*ld,
 98:     where d is the polynomial degree, and X can have ld x d*ld.

100:     Level: advanced

102: .seealso: DSAllocate(), DSGetArray(), DSGetArrayReal(), DSVectors()
103: E*/
104: typedef enum { DS_MAT_A,
105:                DS_MAT_B,
106:                DS_MAT_C,
107:                DS_MAT_T,
108:                DS_MAT_D,
109:                DS_MAT_Q,
110:                DS_MAT_Z,
111:                DS_MAT_X,
112:                DS_MAT_Y,
113:                DS_MAT_U,
114:                DS_MAT_VT,
115:                DS_MAT_W,
116:                DS_MAT_E0,
117:                DS_MAT_E1,
118:                DS_MAT_E2,
119:                DS_MAT_E3,
120:                DS_MAT_E4,
121:                DS_MAT_E5,
122:                DS_MAT_E6,
123:                DS_MAT_E7,
124:                DS_MAT_E8,
125:                DS_MAT_E9,
126:                DS_NUM_MAT } DSMatType;

128: /* Convenience for indexing extra matrices */
129: PETSC_EXTERN DSMatType DSMatExtra[];
130: #define DS_NUM_EXTRA  10

132: PETSC_EXTERN PetscErrorCode DSCreate(MPI_Comm,DS*);
133: PETSC_EXTERN PetscErrorCode DSSetType(DS,DSType);
134: PETSC_EXTERN PetscErrorCode DSGetType(DS,DSType*);
135: PETSC_EXTERN PetscErrorCode DSSetOptionsPrefix(DS,const char *);
136: PETSC_EXTERN PetscErrorCode DSAppendOptionsPrefix(DS,const char *);
137: PETSC_EXTERN PetscErrorCode DSGetOptionsPrefix(DS,const char *[]);
138: PETSC_EXTERN PetscErrorCode DSSetFromOptions(DS);
139: PETSC_EXTERN PetscErrorCode DSView(DS,PetscViewer);
140: PETSC_EXTERN PetscErrorCode DSViewMat(DS,PetscViewer,DSMatType);
141: PETSC_EXTERN PetscErrorCode DSDestroy(DS*);
142: PETSC_EXTERN PetscErrorCode DSReset(DS);

144: PETSC_EXTERN PetscErrorCode DSAllocate(DS,PetscInt);
145: PETSC_EXTERN PetscErrorCode DSGetLeadingDimension(DS,PetscInt*);
146: PETSC_EXTERN PetscErrorCode DSSetState(DS,DSStateType);
147: PETSC_EXTERN PetscErrorCode DSGetState(DS,DSStateType*);
148: PETSC_EXTERN PetscErrorCode DSSetDimensions(DS,PetscInt,PetscInt,PetscInt,PetscInt);
149: PETSC_EXTERN PetscErrorCode DSGetDimensions(DS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
150: PETSC_EXTERN PetscErrorCode DSSetBlockSize(DS,PetscInt);
151: PETSC_EXTERN PetscErrorCode DSGetBlockSize(DS,PetscInt*);
152: PETSC_EXTERN PetscErrorCode DSTruncate(DS,PetscInt);
153: PETSC_EXTERN PetscErrorCode DSSetIdentity(DS,DSMatType);
154: PETSC_EXTERN PetscErrorCode DSSetMethod(DS,PetscInt);
155: PETSC_EXTERN PetscErrorCode DSGetMethod(DS,PetscInt*);
156: PETSC_EXTERN PetscErrorCode DSSetCompact(DS,PetscBool);
157: PETSC_EXTERN PetscErrorCode DSGetCompact(DS,PetscBool*);
158: PETSC_EXTERN PetscErrorCode DSSetExtraRow(DS,PetscBool);
159: PETSC_EXTERN PetscErrorCode DSGetExtraRow(DS,PetscBool*);
160: PETSC_EXTERN PetscErrorCode DSSetRefined(DS,PetscBool);
161: PETSC_EXTERN PetscErrorCode DSGetRefined(DS,PetscBool*);
162: PETSC_EXTERN PetscErrorCode DSGetMat(DS,DSMatType,Mat*);
163: PETSC_EXTERN PetscErrorCode DSRestoreMat(DS,DSMatType,Mat*);
164: PETSC_EXTERN PetscErrorCode DSGetArray(DS,DSMatType,PetscScalar*[]);
165: PETSC_EXTERN PetscErrorCode DSRestoreArray(DS,DSMatType,PetscScalar*[]);
166: PETSC_EXTERN PetscErrorCode DSGetArrayReal(DS,DSMatType,PetscReal*[]);
167: PETSC_EXTERN PetscErrorCode DSRestoreArrayReal(DS,DSMatType,PetscReal*[]);
168: PETSC_EXTERN PetscErrorCode DSVectors(DS,DSMatType,PetscInt*,PetscReal*);
169: PETSC_EXTERN PetscErrorCode DSSolve(DS,PetscScalar*,PetscScalar*);
170: PETSC_EXTERN PetscErrorCode DSSort(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
171: PETSC_EXTERN PetscErrorCode DSCopyMat(DS,DSMatType,PetscInt,PetscInt,Mat,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
172: PETSC_EXTERN PetscErrorCode DSSetSlepcSC(DS,SlepcSC);
173: PETSC_EXTERN PetscErrorCode DSGetSlepcSC(DS,SlepcSC*);
174: PETSC_EXTERN PetscErrorCode DSUpdateExtraRow(DS);
175: PETSC_EXTERN PetscErrorCode DSCond(DS,PetscReal*);
176: PETSC_EXTERN PetscErrorCode DSTranslateHarmonic(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
177: PETSC_EXTERN PetscErrorCode DSTranslateRKS(DS,PetscScalar);
178: PETSC_EXTERN PetscErrorCode DSNormalize(DS,DSMatType,PetscInt);
179: PETSC_EXTERN PetscErrorCode DSOrthogonalize(DS,DSMatType,PetscInt,PetscInt*);
180: PETSC_EXTERN PetscErrorCode DSPseudoOrthogonalize(DS,DSMatType,PetscInt,PetscReal*,PetscInt*,PetscReal*);

182: /* --------- options specific to particular solvers -------- */

184: PETSC_EXTERN PetscErrorCode DSPEPSetDegree(DS,PetscInt);
185: PETSC_EXTERN PetscErrorCode DSPEPGetDegree(DS,PetscInt*);

187: PETSC_EXTERN PetscErrorCode DSNEPSetFN(DS,PetscInt,FN*);
188: PETSC_EXTERN PetscErrorCode DSNEPGetFN(DS,PetscInt,FN*);
189: PETSC_EXTERN PetscErrorCode DSNEPGetNumFN(DS,PetscInt*);

191: PETSC_EXTERN PetscFunctionList DSList;
192: PETSC_EXTERN PetscErrorCode DSRegister(const char[],PetscErrorCode(*)(DS));

194: #endif