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