Actual source code: slepcds.h

slepc-3.17.2 2022-08-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: */

 14: #if !defined(SLEPCDS_H)
 15: #define SLEPCDS_H

 17: #include <slepcsc.h>
 18: #include <slepcfn.h>
 19: #include <slepcrg.h>

 21: #define DS_MAX_SOLVE 6

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

 29:     Level: beginner

 31: .seealso:  DSCreate()
 32: S*/
 33: typedef struct _p_DS* DS;

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

 39:     Level: advanced

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

 55: /* Logging support */
 56: SLEPC_EXTERN PetscClassId DS_CLASSID;

 58: /*E
 59:     DSStateType - Indicates in which state the direct solver is

 61:     Level: advanced

 63: .seealso: DSSetState()
 64: E*/
 65: typedef enum { DS_STATE_RAW,
 66:                DS_STATE_INTERMEDIATE,
 67:                DS_STATE_CONDENSED,
 68:                DS_STATE_TRUNCATED } DSStateType;
 69: SLEPC_EXTERN const char *DSStateTypes[];

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

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

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

 94:     In DSPEP problems, matrices A, B, W can have space for d*ld x d*ld,
 95:     where d is the polynomial degree, and X can have ld x d*ld.
 96:     Also DSNEP has exceptions. Check the manual page of each DS type
 97:     for details.

 99:     Level: advanced

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

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

131: /*E
132:     DSParallelType - Indicates the parallel mode that the direct solver will use

134:     Level: advanced

136: .seealso: DSSetParallel()
137: E*/
138: typedef enum { DS_PARALLEL_REDUNDANT,
139:                DS_PARALLEL_SYNCHRONIZED,
140:                DS_PARALLEL_DISTRIBUTED } DSParallelType;
141: SLEPC_EXTERN const char *DSParallelTypes[];

143: SLEPC_EXTERN PetscErrorCode DSCreate(MPI_Comm,DS*);
144: SLEPC_EXTERN PetscErrorCode DSSetType(DS,DSType);
145: SLEPC_EXTERN PetscErrorCode DSGetType(DS,DSType*);
146: SLEPC_EXTERN PetscErrorCode DSSetOptionsPrefix(DS,const char *);
147: SLEPC_EXTERN PetscErrorCode DSAppendOptionsPrefix(DS,const char *);
148: SLEPC_EXTERN PetscErrorCode DSGetOptionsPrefix(DS,const char *[]);
149: SLEPC_EXTERN PetscErrorCode DSSetFromOptions(DS);
150: SLEPC_EXTERN PetscErrorCode DSView(DS,PetscViewer);
151: SLEPC_EXTERN PetscErrorCode DSViewFromOptions(DS,PetscObject,const char[]);
152: SLEPC_EXTERN PetscErrorCode DSViewMat(DS,PetscViewer,DSMatType);
153: SLEPC_EXTERN PetscErrorCode DSDestroy(DS*);
154: SLEPC_EXTERN PetscErrorCode DSReset(DS);
155: SLEPC_EXTERN PetscErrorCode DSDuplicate(DS,DS*);

157: SLEPC_EXTERN PetscErrorCode DSAllocate(DS,PetscInt);
158: SLEPC_EXTERN PetscErrorCode DSGetLeadingDimension(DS,PetscInt*);
159: SLEPC_EXTERN PetscErrorCode DSSetState(DS,DSStateType);
160: SLEPC_EXTERN PetscErrorCode DSGetState(DS,DSStateType*);
161: SLEPC_EXTERN PetscErrorCode DSSetDimensions(DS,PetscInt,PetscInt,PetscInt);
162: SLEPC_EXTERN PetscErrorCode DSGetDimensions(DS,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
163: SLEPC_EXTERN PetscErrorCode DSSetBlockSize(DS,PetscInt);
164: SLEPC_EXTERN PetscErrorCode DSGetBlockSize(DS,PetscInt*);
165: SLEPC_EXTERN PetscErrorCode DSGetTruncateSize(DS,PetscInt,PetscInt,PetscInt*);
166: SLEPC_EXTERN PetscErrorCode DSTruncate(DS,PetscInt,PetscBool);
167: SLEPC_EXTERN PetscErrorCode DSSetIdentity(DS,DSMatType);
168: SLEPC_EXTERN PetscErrorCode DSSetMethod(DS,PetscInt);
169: SLEPC_EXTERN PetscErrorCode DSGetMethod(DS,PetscInt*);
170: SLEPC_EXTERN PetscErrorCode DSSetParallel(DS,DSParallelType);
171: SLEPC_EXTERN PetscErrorCode DSGetParallel(DS,DSParallelType*);
172: SLEPC_EXTERN PetscErrorCode DSSetCompact(DS,PetscBool);
173: SLEPC_EXTERN PetscErrorCode DSGetCompact(DS,PetscBool*);
174: SLEPC_EXTERN PetscErrorCode DSSetExtraRow(DS,PetscBool);
175: SLEPC_EXTERN PetscErrorCode DSGetExtraRow(DS,PetscBool*);
176: SLEPC_EXTERN PetscErrorCode DSSetRefined(DS,PetscBool);
177: SLEPC_EXTERN PetscErrorCode DSGetRefined(DS,PetscBool*);
178: SLEPC_EXTERN PetscErrorCode DSGetMat(DS,DSMatType,Mat*);
179: SLEPC_EXTERN PetscErrorCode DSRestoreMat(DS,DSMatType,Mat*);
180: SLEPC_EXTERN PetscErrorCode DSGetArray(DS,DSMatType,PetscScalar*[]);
181: SLEPC_EXTERN PetscErrorCode DSRestoreArray(DS,DSMatType,PetscScalar*[]);
182: SLEPC_EXTERN PetscErrorCode DSGetArrayReal(DS,DSMatType,PetscReal*[]);
183: SLEPC_EXTERN PetscErrorCode DSRestoreArrayReal(DS,DSMatType,PetscReal*[]);
184: SLEPC_EXTERN PetscErrorCode DSVectors(DS,DSMatType,PetscInt*,PetscReal*);
185: SLEPC_EXTERN PetscErrorCode DSSolve(DS,PetscScalar*,PetscScalar*);
186: SLEPC_EXTERN PetscErrorCode DSSort(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
187: SLEPC_EXTERN PetscErrorCode DSSortWithPermutation(DS,PetscInt*,PetscScalar*,PetscScalar*);
188: SLEPC_EXTERN PetscErrorCode DSSynchronize(DS,PetscScalar*,PetscScalar*);
189: SLEPC_EXTERN PetscErrorCode DSCopyMat(DS,DSMatType,PetscInt,PetscInt,Mat,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
190: SLEPC_EXTERN PetscErrorCode DSMatGetSize(DS,DSMatType,PetscInt*,PetscInt*);
191: SLEPC_EXTERN PetscErrorCode DSMatIsHermitian(DS,DSMatType,PetscBool*);
192: SLEPC_EXTERN PetscErrorCode DSSetSlepcSC(DS,SlepcSC);
193: SLEPC_EXTERN PetscErrorCode DSGetSlepcSC(DS,SlepcSC*);
194: SLEPC_EXTERN PetscErrorCode DSUpdateExtraRow(DS);
195: SLEPC_EXTERN PetscErrorCode DSCond(DS,PetscReal*);
196: SLEPC_EXTERN PetscErrorCode DSTranslateHarmonic(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
197: SLEPC_EXTERN PetscErrorCode DSTranslateRKS(DS,PetscScalar);
198: SLEPC_EXTERN PetscErrorCode DSOrthogonalize(DS,DSMatType,PetscInt,PetscInt*);
199: SLEPC_EXTERN PetscErrorCode DSPseudoOrthogonalize(DS,DSMatType,PetscInt,PetscReal*,PetscInt*,PetscReal*);

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

203: SLEPC_EXTERN PetscErrorCode DSSVDSetDimensions(DS,PetscInt);
204: SLEPC_EXTERN PetscErrorCode DSSVDGetDimensions(DS,PetscInt*);
205: SLEPC_EXTERN PetscErrorCode DSGSVDSetDimensions(DS,PetscInt,PetscInt);
206: SLEPC_EXTERN PetscErrorCode DSGSVDGetDimensions(DS,PetscInt*,PetscInt*);

208: SLEPC_EXTERN PetscErrorCode DSPEPSetDegree(DS,PetscInt);
209: SLEPC_EXTERN PetscErrorCode DSPEPGetDegree(DS,PetscInt*);
210: SLEPC_EXTERN PetscErrorCode DSPEPSetCoefficients(DS,PetscReal*);
211: SLEPC_EXTERN PetscErrorCode DSPEPGetCoefficients(DS,PetscReal**);

213: SLEPC_EXTERN PetscErrorCode DSNEPSetFN(DS,PetscInt,FN*);
214: SLEPC_EXTERN PetscErrorCode DSNEPGetFN(DS,PetscInt,FN*);
215: SLEPC_EXTERN PetscErrorCode DSNEPGetNumFN(DS,PetscInt*);
216: SLEPC_EXTERN PetscErrorCode DSNEPSetMinimality(DS,PetscInt);
217: SLEPC_EXTERN PetscErrorCode DSNEPGetMinimality(DS,PetscInt*);
218: SLEPC_EXTERN PetscErrorCode DSNEPSetRefine(DS,PetscReal,PetscInt);
219: SLEPC_EXTERN PetscErrorCode DSNEPGetRefine(DS,PetscReal*,PetscInt*);
220: SLEPC_EXTERN PetscErrorCode DSNEPSetIntegrationPoints(DS,PetscInt);
221: SLEPC_EXTERN PetscErrorCode DSNEPGetIntegrationPoints(DS,PetscInt*);
222: SLEPC_EXTERN PetscErrorCode DSNEPSetSamplingSize(DS,PetscInt);
223: SLEPC_EXTERN PetscErrorCode DSNEPGetSamplingSize(DS,PetscInt*);
224: SLEPC_EXTERN PetscErrorCode DSNEPSetRG(DS,RG);
225: SLEPC_EXTERN PetscErrorCode DSNEPGetRG(DS,RG*);
226: SLEPC_EXTERN PetscErrorCode DSNEPSetComputeMatrixFunction(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*);
227: SLEPC_EXTERN PetscErrorCode DSNEPGetComputeMatrixFunction(DS,PetscErrorCode (**)(DS,PetscScalar,PetscBool,DSMatType,void*),void**);

229: SLEPC_EXTERN PetscFunctionList DSList;
230: SLEPC_EXTERN PetscErrorCode DSRegister(const char[],PetscErrorCode(*)(DS));

232: #endif