Actual source code: slepcds.h

slepc-3.22.1 2024-10-28
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: #pragma once

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

 20: /* SUBMANSEC = DS */

 22: #define DS_MAX_SOLVE 6

 24: SLEPC_EXTERN PetscErrorCode DSInitializePackage(void);
 25: SLEPC_EXTERN PetscErrorCode DSFinalizePackage(void);

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

 32:     Level: beginner

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

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

 42:     Level: advanced

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

 59: /* Logging support */
 60: SLEPC_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;
 73: SLEPC_EXTERN const char *DSStateTypes[];

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

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

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

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

103:     Level: advanced

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

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

135: /*E
136:     DSParallelType - Indicates the parallel mode that the direct solver will use

138:     Level: advanced

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

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

161: SLEPC_EXTERN PetscErrorCode DSAllocate(DS,PetscInt);
162: SLEPC_EXTERN PetscErrorCode DSGetLeadingDimension(DS,PetscInt*);
163: SLEPC_EXTERN PetscErrorCode DSSetState(DS,DSStateType);
164: SLEPC_EXTERN PetscErrorCode DSGetState(DS,DSStateType*);
165: SLEPC_EXTERN PetscErrorCode DSSetDimensions(DS,PetscInt,PetscInt,PetscInt);
166: SLEPC_EXTERN PetscErrorCode DSGetDimensions(DS,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
167: SLEPC_EXTERN PetscErrorCode DSSetBlockSize(DS,PetscInt);
168: SLEPC_EXTERN PetscErrorCode DSGetBlockSize(DS,PetscInt*);
169: SLEPC_EXTERN PetscErrorCode DSGetTruncateSize(DS,PetscInt,PetscInt,PetscInt*);
170: SLEPC_EXTERN PetscErrorCode DSTruncate(DS,PetscInt,PetscBool);
171: SLEPC_EXTERN PetscErrorCode DSSetIdentity(DS,DSMatType);
172: SLEPC_EXTERN PetscErrorCode DSSetMethod(DS,PetscInt);
173: SLEPC_EXTERN PetscErrorCode DSGetMethod(DS,PetscInt*);
174: SLEPC_EXTERN PetscErrorCode DSSetParallel(DS,DSParallelType);
175: SLEPC_EXTERN PetscErrorCode DSGetParallel(DS,DSParallelType*);
176: SLEPC_EXTERN PetscErrorCode DSSetCompact(DS,PetscBool);
177: SLEPC_EXTERN PetscErrorCode DSGetCompact(DS,PetscBool*);
178: SLEPC_EXTERN PetscErrorCode DSSetExtraRow(DS,PetscBool);
179: SLEPC_EXTERN PetscErrorCode DSGetExtraRow(DS,PetscBool*);
180: SLEPC_EXTERN PetscErrorCode DSSetRefined(DS,PetscBool);
181: SLEPC_EXTERN PetscErrorCode DSGetRefined(DS,PetscBool*);
182: SLEPC_EXTERN PetscErrorCode DSGetMat(DS,DSMatType,Mat*);
183: SLEPC_EXTERN PetscErrorCode DSRestoreMat(DS,DSMatType,Mat*);
184: SLEPC_EXTERN PetscErrorCode DSGetMatAndColumn(DS,DSMatType,PetscInt,Mat*,Vec*);
185: SLEPC_EXTERN PetscErrorCode DSRestoreMatAndColumn(DS,DSMatType,PetscInt,Mat*,Vec*);
186: SLEPC_EXTERN PetscErrorCode DSGetArray(DS,DSMatType,PetscScalar*[]);
187: SLEPC_EXTERN PetscErrorCode DSRestoreArray(DS,DSMatType,PetscScalar*[]);
188: SLEPC_EXTERN PetscErrorCode DSGetArrayReal(DS,DSMatType,PetscReal*[]);
189: SLEPC_EXTERN PetscErrorCode DSRestoreArrayReal(DS,DSMatType,PetscReal*[]);
190: SLEPC_EXTERN PetscErrorCode DSVectors(DS,DSMatType,PetscInt*,PetscReal*);
191: SLEPC_EXTERN PetscErrorCode DSSolve(DS,PetscScalar*,PetscScalar*);
192: SLEPC_EXTERN PetscErrorCode DSSort(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
193: SLEPC_EXTERN PetscErrorCode DSSortWithPermutation(DS,PetscInt*,PetscScalar*,PetscScalar*);
194: SLEPC_EXTERN PetscErrorCode DSSynchronize(DS,PetscScalar*,PetscScalar*);
195: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "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)
196: {
197:   Mat M,M0,A0;

199:   PetscFunctionBegin;
200:   PetscCall(DSGetMat(ds,m,&M));
201:   PetscCall(MatDenseGetSubMatrix(M,mr,mr+rows,mc,mc+cols,&M0));
202:   PetscCall(MatDenseGetSubMatrix(A,Ar,Ar+rows,Ac,Ac+cols,&A0));
203:   if (out) PetscCall(MatCopy(M0,A0,SAME_NONZERO_PATTERN));
204:   else PetscCall(MatCopy(A0,M0,SAME_NONZERO_PATTERN));
205:   PetscCall(MatDenseRestoreSubMatrix(M,&M0));
206:   PetscCall(MatDenseRestoreSubMatrix(A,&A0));
207:   PetscCall(DSRestoreMat(ds,m,&M));
208:   PetscFunctionReturn(PETSC_SUCCESS);
209: }
210: SLEPC_EXTERN PetscErrorCode DSMatGetSize(DS,DSMatType,PetscInt*,PetscInt*);
211: SLEPC_EXTERN PetscErrorCode DSMatIsHermitian(DS,DSMatType,PetscBool*);
212: SLEPC_EXTERN PetscErrorCode DSSetSlepcSC(DS,SlepcSC);
213: SLEPC_EXTERN PetscErrorCode DSGetSlepcSC(DS,SlepcSC*);
214: SLEPC_EXTERN PetscErrorCode DSUpdateExtraRow(DS);
215: SLEPC_EXTERN PetscErrorCode DSCond(DS,PetscReal*);
216: SLEPC_EXTERN PetscErrorCode DSTranslateHarmonic(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
217: SLEPC_EXTERN PetscErrorCode DSTranslateRKS(DS,PetscScalar);
218: SLEPC_EXTERN PetscErrorCode DSOrthogonalize(DS,DSMatType,PetscInt,PetscInt*);
219: SLEPC_EXTERN PetscErrorCode DSPseudoOrthogonalize(DS,DSMatType,PetscInt,PetscReal*,PetscInt*,PetscReal*);

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

223: SLEPC_EXTERN PetscErrorCode DSSVDSetDimensions(DS,PetscInt);
224: SLEPC_EXTERN PetscErrorCode DSSVDGetDimensions(DS,PetscInt*);
225: SLEPC_EXTERN PetscErrorCode DSGSVDSetDimensions(DS,PetscInt,PetscInt);
226: SLEPC_EXTERN PetscErrorCode DSGSVDGetDimensions(DS,PetscInt*,PetscInt*);
227: SLEPC_EXTERN PetscErrorCode DSHSVDSetDimensions(DS,PetscInt);
228: SLEPC_EXTERN PetscErrorCode DSHSVDGetDimensions(DS,PetscInt*);
229: SLEPC_EXTERN PetscErrorCode DSHSVDSetReorthogonalize(DS,PetscBool);
230: SLEPC_EXTERN PetscErrorCode DSHSVDGetReorthogonalize(DS,PetscBool*);

232: SLEPC_EXTERN PetscErrorCode DSPEPSetDegree(DS,PetscInt);
233: SLEPC_EXTERN PetscErrorCode DSPEPGetDegree(DS,PetscInt*);
234: SLEPC_EXTERN PetscErrorCode DSPEPSetCoefficients(DS,PetscReal*);
235: SLEPC_EXTERN PetscErrorCode DSPEPGetCoefficients(DS,PetscReal**);

237: SLEPC_EXTERN PetscErrorCode DSNEPSetFN(DS,PetscInt,FN*);
238: SLEPC_EXTERN PetscErrorCode DSNEPGetFN(DS,PetscInt,FN*);
239: SLEPC_EXTERN PetscErrorCode DSNEPGetNumFN(DS,PetscInt*);
240: SLEPC_EXTERN PetscErrorCode DSNEPSetMinimality(DS,PetscInt);
241: SLEPC_EXTERN PetscErrorCode DSNEPGetMinimality(DS,PetscInt*);
242: SLEPC_EXTERN PetscErrorCode DSNEPSetRefine(DS,PetscReal,PetscInt);
243: SLEPC_EXTERN PetscErrorCode DSNEPGetRefine(DS,PetscReal*,PetscInt*);
244: SLEPC_EXTERN PetscErrorCode DSNEPSetIntegrationPoints(DS,PetscInt);
245: SLEPC_EXTERN PetscErrorCode DSNEPGetIntegrationPoints(DS,PetscInt*);
246: SLEPC_EXTERN PetscErrorCode DSNEPSetSamplingSize(DS,PetscInt);
247: SLEPC_EXTERN PetscErrorCode DSNEPGetSamplingSize(DS,PetscInt*);
248: SLEPC_EXTERN PetscErrorCode DSNEPSetRG(DS,RG);
249: SLEPC_EXTERN PetscErrorCode DSNEPGetRG(DS,RG*);

251: /*S
252:   DSNEPMatrixFunctionFn - A prototype of a DSNEP compute matrix function that would be passed to DSNEPSetComputeMatrixFunction()

254:   Calling Sequence:
255: +   ds     - the direct solver object
256: .   lambda - point where T(lambda) or T'(lambda) must be evaluated
257: .   deriv  - if true compute T'(lambda), otherwise compute T(lambda)
258: .   mat    - the DS matrix where the result must be stored
259: -   ctx    - [optional] user-defined context for private data for the
260:              matrix evaluation routine (may be `NULL`)

262:   Level: developer

264: .seealso: DSNEPSetComputeMatrixFunction()
265: S*/
266: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DSNEPMatrixFunctionFn)(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx);

268: SLEPC_EXTERN PetscErrorCode DSNEPSetComputeMatrixFunction(DS,DSNEPMatrixFunctionFn*,void*);
269: SLEPC_EXTERN PetscErrorCode DSNEPGetComputeMatrixFunction(DS,DSNEPMatrixFunctionFn**,void**);

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