Actual source code: dsimpl.h
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: */
11: #if !defined(SLEPCDSIMPL_H)
12: #define SLEPCDSIMPL_H
14: #include <slepcds.h>
15: #include <slepc/private/slepcimpl.h>
17: SLEPC_EXTERN PetscBool DSRegisterAllCalled;
18: SLEPC_EXTERN PetscErrorCode DSRegisterAll(void);
19: SLEPC_EXTERN PetscLogEvent DS_Solve,DS_Vectors,DS_Synchronize,DS_Other;
20: SLEPC_INTERN const char *DSMatName[];
22: typedef struct _DSOps *DSOps;
24: struct _DSOps {
25: PetscErrorCode (*allocate)(DS,PetscInt);
26: PetscErrorCode (*view)(DS,PetscViewer);
27: PetscErrorCode (*vectors)(DS,DSMatType,PetscInt*,PetscReal*);
28: PetscErrorCode (*solve[DS_MAX_SOLVE])(DS,PetscScalar*,PetscScalar*);
29: PetscErrorCode (*sort)(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
30: PetscErrorCode (*sortperm)(DS,PetscInt*,PetscScalar*,PetscScalar*);
31: PetscErrorCode (*truncate)(DS,PetscInt);
32: PetscErrorCode (*update)(DS);
33: PetscErrorCode (*cond)(DS,PetscReal*);
34: PetscErrorCode (*transharm)(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
35: PetscErrorCode (*transrks)(DS,PetscScalar);
36: PetscErrorCode (*destroy)(DS);
37: PetscErrorCode (*matgetsize)(DS,DSMatType,PetscInt*,PetscInt*);
38: PetscErrorCode (*hermitian)(DS,DSMatType,PetscBool*);
39: PetscErrorCode (*synchronize)(DS,PetscScalar*,PetscScalar*);
40: };
42: struct _p_DS {
43: PETSCHEADER(struct _DSOps);
44: /*------------------------- User parameters --------------------------*/
45: DSStateType state; /* the current state */
46: PetscInt method; /* identifies the variant to be used */
47: PetscBool compact; /* whether the matrices are stored in compact form */
48: PetscBool refined; /* get refined vectors instead of regular vectors */
49: PetscBool extrarow; /* assume the matrix dimension is (n+1) x n */
50: PetscInt ld; /* leading dimension */
51: PetscInt l; /* number of locked (inactive) leading columns */
52: PetscInt n; /* current dimension */
53: PetscInt m; /* current column dimension (for SVD only) */
54: PetscInt k; /* intermediate dimension (e.g. position of arrow) */
55: PetscInt t; /* length of decomposition when it was truncated */
56: PetscInt bs; /* block size */
57: SlepcSC sc; /* sorting criterion */
58: DSParallelType pmode; /* parallel mode (redundant or synchronized) */
60: /*----------------- Status variables and working data ----------------*/
61: PetscScalar *mat[DS_NUM_MAT]; /* the matrices */
62: PetscReal *rmat[DS_NUM_MAT]; /* the matrices (real) */
63: Mat omat[DS_NUM_MAT]; /* the matrices (PETSc object) */
64: PetscInt *perm; /* permutation */
65: void *data; /* placeholder for solver-specific stuff */
66: PetscBool scset; /* the sc was provided by the user */
67: PetscScalar *work;
68: PetscReal *rwork;
69: PetscBLASInt *iwork;
70: PetscInt lwork,lrwork,liwork;
71: };
73: /*
74: Macros to test valid DS arguments
75: */
76: #if !defined(PETSC_USE_DEBUG)
78: #define DSCheckAlloc(h,arg) do {} while (0)
79: #define DSCheckSolved(h,arg) do {} while (0)
80: #define DSCheckValidMat(ds,m,arg) do {} while (0)
81: #define DSCheckValidMatReal(ds,m,arg) do {} while (0)
83: #else
85: #define DSCheckAlloc(h,arg) \
86: do { \
87: if (!(h)->ld) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call DSAllocate() first: Parameter #%d",arg); \
88: } while (0)
90: #define DSCheckSolved(h,arg) \
91: do { \
92: if ((h)->state<DS_STATE_CONDENSED) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call DSSolve() first: Parameter #%d",arg); \
93: } while (0)
95: #define DSCheckValidMat(ds,m,arg) \
96: do { \
97: if ((m)>=DS_NUM_MAT) SETERRQ1(PetscObjectComm((PetscObject)(ds)),PETSC_ERR_ARG_WRONG,"Invalid matrix: Parameter #%d",arg); \
98: if (!(ds)->mat[m]) SETERRQ1(PetscObjectComm((PetscObject)(ds)),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS: Parameter #%d",arg); \
99: } while (0)
101: #define DSCheckValidMatReal(ds,m,arg) \
102: do { \
103: if ((m)>=DS_NUM_MAT) SETERRQ1(PetscObjectComm((PetscObject)(ds)),PETSC_ERR_ARG_WRONG,"Invalid matrix: Parameter #%d",arg); \
104: if (!(ds)->rmat[m]) SETERRQ1(PetscObjectComm((PetscObject)(ds)),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS: Parameter #%d",arg); \
105: } while (0)
107: #endif
109: SLEPC_INTERN PetscErrorCode DSAllocateMatrix_Private(DS,DSMatType,PetscBool);
110: PETSC_STATIC_INLINE PetscErrorCode DSAllocateMat_Private(DS ds,DSMatType m) {return DSAllocateMatrix_Private(ds,m,PETSC_FALSE);}
111: PETSC_STATIC_INLINE PetscErrorCode DSAllocateMatReal_Private(DS ds,DSMatType m) {return DSAllocateMatrix_Private(ds,m,PETSC_TRUE);}
112: SLEPC_INTERN PetscErrorCode DSAllocateWork_Private(DS,PetscInt,PetscInt,PetscInt);
113: SLEPC_INTERN PetscErrorCode DSSortEigenvalues_Private(DS,PetscScalar*,PetscScalar*,PetscInt*,PetscBool);
114: SLEPC_INTERN PetscErrorCode DSSortEigenvaluesReal_Private(DS,PetscReal*,PetscInt*);
115: SLEPC_INTERN PetscErrorCode DSPermuteColumns_Private(DS,PetscInt,PetscInt,DSMatType,PetscInt*);
116: SLEPC_INTERN PetscErrorCode DSPermuteRows_Private(DS,PetscInt,PetscInt,DSMatType,PetscInt*);
117: SLEPC_INTERN PetscErrorCode DSPermuteBoth_Private(DS,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
118: SLEPC_INTERN PetscErrorCode DSCopyMatrix_Private(DS,DSMatType,DSMatType);
120: SLEPC_INTERN PetscErrorCode DSGHIEPOrthogEigenv(DS,DSMatType,PetscScalar*,PetscScalar*,PetscBool);
121: SLEPC_INTERN PetscErrorCode DSGHIEPComplexEigs(DS,PetscInt,PetscInt,PetscScalar*,PetscScalar*);
122: SLEPC_INTERN PetscErrorCode DSGHIEPInverseIteration(DS,PetscScalar*,PetscScalar*);
123: SLEPC_INTERN PetscErrorCode DSIntermediate_GHIEP(DS);
124: SLEPC_INTERN PetscErrorCode DSSwitchFormat_GHIEP(DS,PetscBool);
125: SLEPC_INTERN PetscErrorCode DSGHIEPRealBlocks(DS);
126: SLEPC_INTERN PetscErrorCode DSSolve_GHIEP_HZ(DS,PetscScalar*,PetscScalar*);
128: SLEPC_INTERN PetscErrorCode BDC_dibtdc_(const char*,PetscBLASInt,PetscBLASInt,PetscBLASInt*,PetscReal*,PetscBLASInt,PetscBLASInt,PetscReal*,PetscBLASInt*,PetscBLASInt,PetscBLASInt,PetscReal,PetscReal*,PetscReal*,PetscBLASInt,PetscReal*,PetscBLASInt,PetscBLASInt*,PetscBLASInt,PetscBLASInt*,PetscBLASInt);
129: SLEPC_INTERN PetscErrorCode BDC_dlaed3m_(const char*,const char*,PetscBLASInt,PetscBLASInt,PetscBLASInt,PetscReal*,PetscReal*,PetscBLASInt,PetscReal,PetscReal*,PetscReal*,PetscBLASInt*,PetscBLASInt*,PetscReal*,PetscReal*,PetscBLASInt*,PetscBLASInt,PetscBLASInt);
130: SLEPC_INTERN PetscErrorCode BDC_dmerg2_(const char*,PetscBLASInt,PetscBLASInt,PetscReal*,PetscReal*,PetscBLASInt,PetscBLASInt*,PetscReal*,PetscReal*,PetscBLASInt,PetscReal*,PetscBLASInt,PetscBLASInt,PetscReal*,PetscBLASInt,PetscBLASInt*,PetscReal,PetscBLASInt*,PetscBLASInt);
131: SLEPC_INTERN PetscErrorCode BDC_dsbtdc_(const char*,const char*,PetscBLASInt,PetscBLASInt,PetscBLASInt*,PetscReal*,PetscBLASInt,PetscBLASInt,PetscReal*,PetscBLASInt,PetscBLASInt,PetscReal,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBLASInt,PetscReal*,PetscBLASInt,PetscBLASInt*,PetscBLASInt,PetscReal*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt,PetscBLASInt);
132: SLEPC_INTERN PetscErrorCode BDC_dsrtdf_(PetscBLASInt*,PetscBLASInt,PetscBLASInt,PetscReal*,PetscReal*,PetscBLASInt,PetscBLASInt*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscReal,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
134: #endif