Actual source code: slepcds.h

slepc-3.18.1 2022-11-02
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: /* SUBMANSEC = DS */

 23: #define DS_MAX_SOLVE 6

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

 31:     Level: beginner

 33: .seealso:  DSCreate()
 34: S*/
 35: typedef struct _p_DS* DS;

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

 41:     Level: advanced

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

 57: /* Logging support */
 58: SLEPC_EXTERN PetscClassId DS_CLASSID;

 60: /*E
 61:     DSStateType - Indicates in which state the direct solver is

 63:     Level: advanced

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

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

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

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

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

101:     Level: advanced

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

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

133: /*E
134:     DSParallelType - Indicates the parallel mode that the direct solver will use

136:     Level: advanced

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

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

159: SLEPC_EXTERN PetscErrorCode DSAllocate(DS,PetscInt);
160: SLEPC_EXTERN PetscErrorCode DSGetLeadingDimension(DS,PetscInt*);
161: SLEPC_EXTERN PetscErrorCode DSSetState(DS,DSStateType);
162: SLEPC_EXTERN PetscErrorCode DSGetState(DS,DSStateType*);
163: SLEPC_EXTERN PetscErrorCode DSSetDimensions(DS,PetscInt,PetscInt,PetscInt);
164: SLEPC_EXTERN PetscErrorCode DSGetDimensions(DS,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
165: SLEPC_EXTERN PetscErrorCode DSSetBlockSize(DS,PetscInt);
166: SLEPC_EXTERN PetscErrorCode DSGetBlockSize(DS,PetscInt*);
167: SLEPC_EXTERN PetscErrorCode DSGetTruncateSize(DS,PetscInt,PetscInt,PetscInt*);
168: SLEPC_EXTERN PetscErrorCode DSTruncate(DS,PetscInt,PetscBool);
169: SLEPC_EXTERN PetscErrorCode DSSetIdentity(DS,DSMatType);
170: SLEPC_EXTERN PetscErrorCode DSSetMethod(DS,PetscInt);
171: SLEPC_EXTERN PetscErrorCode DSGetMethod(DS,PetscInt*);
172: SLEPC_EXTERN PetscErrorCode DSSetParallel(DS,DSParallelType);
173: SLEPC_EXTERN PetscErrorCode DSGetParallel(DS,DSParallelType*);
174: SLEPC_EXTERN PetscErrorCode DSSetCompact(DS,PetscBool);
175: SLEPC_EXTERN PetscErrorCode DSGetCompact(DS,PetscBool*);
176: SLEPC_EXTERN PetscErrorCode DSSetExtraRow(DS,PetscBool);
177: SLEPC_EXTERN PetscErrorCode DSGetExtraRow(DS,PetscBool*);
178: SLEPC_EXTERN PetscErrorCode DSSetRefined(DS,PetscBool);
179: SLEPC_EXTERN PetscErrorCode DSGetRefined(DS,PetscBool*);
180: SLEPC_EXTERN PetscErrorCode DSGetMat(DS,DSMatType,Mat*);
181: SLEPC_EXTERN PetscErrorCode DSRestoreMat(DS,DSMatType,Mat*);
182: SLEPC_EXTERN PetscErrorCode DSGetMatAndColumn(DS,DSMatType,PetscInt,Mat*,Vec*);
183: SLEPC_EXTERN PetscErrorCode DSRestoreMatAndColumn(DS,DSMatType,PetscInt,Mat*,Vec*);
184: SLEPC_EXTERN PetscErrorCode DSGetArray(DS,DSMatType,PetscScalar*[]);
185: SLEPC_EXTERN PetscErrorCode DSRestoreArray(DS,DSMatType,PetscScalar*[]);
186: SLEPC_EXTERN PetscErrorCode DSGetArrayReal(DS,DSMatType,PetscReal*[]);
187: SLEPC_EXTERN PetscErrorCode DSRestoreArrayReal(DS,DSMatType,PetscReal*[]);
188: SLEPC_EXTERN PetscErrorCode DSVectors(DS,DSMatType,PetscInt*,PetscReal*);
189: SLEPC_EXTERN PetscErrorCode DSSolve(DS,PetscScalar*,PetscScalar*);
190: SLEPC_EXTERN PetscErrorCode DSSort(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
191: SLEPC_EXTERN PetscErrorCode DSSortWithPermutation(DS,PetscInt*,PetscScalar*,PetscScalar*);
192: SLEPC_EXTERN PetscErrorCode DSSynchronize(DS,PetscScalar*,PetscScalar*);
193: PETSC_DEPRECATED_FUNCTION("Use DSGetMat()+MatDenseGetSubMatrix()+MatCopy()") static inline PetscErrorCode DSCopyMat(DS ds,DSMatType m,PetscInt mr,PetscInt mc,Mat A,PetscInt Ar,PetscInt Ac,PetscInt rows,PetscInt cols,PetscBool out)
194: {
195:   Mat M,M0,A0;
196:   DSGetMat(ds,m,&M);
197:   MatDenseGetSubMatrix(M,mr,mr+rows,mc,mc+cols,&M0);
198:   MatDenseGetSubMatrix(A,Ar,Ar+rows,Ac,Ac+cols,&A0);
199:   if (out) MatCopy(M0,A0,SAME_NONZERO_PATTERN);
200:   else MatCopy(A0,M0,SAME_NONZERO_PATTERN);
201:   MatDenseRestoreSubMatrix(M,&M0);
202:   MatDenseRestoreSubMatrix(A,&A0);
203:   DSRestoreMat(ds,m,&M);
204:   return 0;
205: }
206: SLEPC_EXTERN PetscErrorCode DSMatGetSize(DS,DSMatType,PetscInt*,PetscInt*);
207: SLEPC_EXTERN PetscErrorCode DSMatIsHermitian(DS,DSMatType,PetscBool*);
208: SLEPC_EXTERN PetscErrorCode DSSetSlepcSC(DS,SlepcSC);
209: SLEPC_EXTERN PetscErrorCode DSGetSlepcSC(DS,SlepcSC*);
210: SLEPC_EXTERN PetscErrorCode DSUpdateExtraRow(DS);
211: SLEPC_EXTERN PetscErrorCode DSCond(DS,PetscReal*);
212: SLEPC_EXTERN PetscErrorCode DSTranslateHarmonic(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
213: SLEPC_EXTERN PetscErrorCode DSTranslateRKS(DS,PetscScalar);
214: SLEPC_EXTERN PetscErrorCode DSOrthogonalize(DS,DSMatType,PetscInt,PetscInt*);
215: SLEPC_EXTERN PetscErrorCode DSPseudoOrthogonalize(DS,DSMatType,PetscInt,PetscReal*,PetscInt*,PetscReal*);

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

219: SLEPC_EXTERN PetscErrorCode DSSVDSetDimensions(DS,PetscInt);
220: SLEPC_EXTERN PetscErrorCode DSSVDGetDimensions(DS,PetscInt*);
221: SLEPC_EXTERN PetscErrorCode DSGSVDSetDimensions(DS,PetscInt,PetscInt);
222: SLEPC_EXTERN PetscErrorCode DSGSVDGetDimensions(DS,PetscInt*,PetscInt*);

224: SLEPC_EXTERN PetscErrorCode DSPEPSetDegree(DS,PetscInt);
225: SLEPC_EXTERN PetscErrorCode DSPEPGetDegree(DS,PetscInt*);
226: SLEPC_EXTERN PetscErrorCode DSPEPSetCoefficients(DS,PetscReal*);
227: SLEPC_EXTERN PetscErrorCode DSPEPGetCoefficients(DS,PetscReal**);

229: SLEPC_EXTERN PetscErrorCode DSNEPSetFN(DS,PetscInt,FN*);
230: SLEPC_EXTERN PetscErrorCode DSNEPGetFN(DS,PetscInt,FN*);
231: SLEPC_EXTERN PetscErrorCode DSNEPGetNumFN(DS,PetscInt*);
232: SLEPC_EXTERN PetscErrorCode DSNEPSetMinimality(DS,PetscInt);
233: SLEPC_EXTERN PetscErrorCode DSNEPGetMinimality(DS,PetscInt*);
234: SLEPC_EXTERN PetscErrorCode DSNEPSetRefine(DS,PetscReal,PetscInt);
235: SLEPC_EXTERN PetscErrorCode DSNEPGetRefine(DS,PetscReal*,PetscInt*);
236: SLEPC_EXTERN PetscErrorCode DSNEPSetIntegrationPoints(DS,PetscInt);
237: SLEPC_EXTERN PetscErrorCode DSNEPGetIntegrationPoints(DS,PetscInt*);
238: SLEPC_EXTERN PetscErrorCode DSNEPSetSamplingSize(DS,PetscInt);
239: SLEPC_EXTERN PetscErrorCode DSNEPGetSamplingSize(DS,PetscInt*);
240: SLEPC_EXTERN PetscErrorCode DSNEPSetRG(DS,RG);
241: SLEPC_EXTERN PetscErrorCode DSNEPGetRG(DS,RG*);
242: SLEPC_EXTERN PetscErrorCode DSNEPSetComputeMatrixFunction(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*);
243: SLEPC_EXTERN PetscErrorCode DSNEPGetComputeMatrixFunction(DS,PetscErrorCode (**)(DS,PetscScalar,PetscBool,DSMatType,void*),void**);

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

248: #endif