Actual source code: slepcds.h

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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 the direct solver object in SLEPc
 12: */

 16:  #include <slepcsc.h>
 17:  #include <slepcfn.h>

 19: #define DS_MAX_SOLVE 6

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

 27:     Level: beginner

 29: .seealso:  DSCreate()
 30: S*/
 31: typedef struct _p_DS* DS;

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

 37:     Level: advanced

 39: .seealso: DSSetType(), DS
 40: J*/
 41: typedef const char* DSType;
 42: #define DSHEP    "hep"
 43: #define DSNHEP   "nhep"
 44: #define DSGHEP   "ghep"
 45: #define DSGHIEP  "ghiep"
 46: #define DSGNHEP  "gnhep"
 47: #define DSSVD    "svd"
 48: #define DSPEP    "pep"
 49: #define DSNEP    "nep"

 51: /* Logging support */
 52: PETSC_EXTERN PetscClassId DS_CLASSID;

 54: /*E
 55:     DSStateType - Indicates in which state the direct solver is

 57:     Level: advanced

 59: .seealso: DSSetState()
 60: E*/
 61: typedef enum { DS_STATE_RAW,
 62:                DS_STATE_INTERMEDIATE,
 63:                DS_STATE_CONDENSED,
 64:                DS_STATE_TRUNCATED } DSStateType;
 65: PETSC_EXTERN const char *DSStateTypes[];

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

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

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

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

 93:     Level: advanced

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

121: /* Convenience for indexing extra matrices */
122: PETSC_EXTERN DSMatType DSMatExtra[];
123: #define DS_NUM_EXTRA  10

125: /*E
126:     DSParallelType - Indicates the parallel mode that the direct solver will use

128:     Level: advanced

130: .seealso: DSSetParallel()
131: E*/
132: typedef enum { DS_PARALLEL_REDUNDANT,
133:                DS_PARALLEL_SYNCHRONIZED } DSParallelType;
134: PETSC_EXTERN const char *DSParallelTypes[];

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

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

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

192: PETSC_EXTERN PetscErrorCode DSPEPSetDegree(DS,PetscInt);
193: PETSC_EXTERN PetscErrorCode DSPEPGetDegree(DS,PetscInt*);

195: PETSC_EXTERN PetscErrorCode DSNEPSetFN(DS,PetscInt,FN*);
196: PETSC_EXTERN PetscErrorCode DSNEPGetFN(DS,PetscInt,FN*);
197: PETSC_EXTERN PetscErrorCode DSNEPGetNumFN(DS,PetscInt*);

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

202: #endif