Actual source code: dsimpl.h

slepc-3.18.0 2022-10-01
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: */

 11: #if !defined(SLEPCDSIMPL_H)
 12: #define SLEPCDSIMPL_H

 14: #include <slepcds.h>
 15: #include <slepc/private/slepcimpl.h>

 17: /* SUBMANSEC = DS */

 19: SLEPC_EXTERN PetscBool DSRegisterAllCalled;
 20: SLEPC_EXTERN PetscErrorCode DSRegisterAll(void);
 21: SLEPC_EXTERN PetscLogEvent DS_Solve,DS_Vectors,DS_Synchronize,DS_Other;
 22: SLEPC_INTERN const char *DSMatName[];

 24: typedef struct _DSOps *DSOps;

 26: struct _DSOps {
 27:   PetscErrorCode (*allocate)(DS,PetscInt);
 28:   PetscErrorCode (*setfromoptions)(DS,PetscOptionItems*);
 29:   PetscErrorCode (*view)(DS,PetscViewer);
 30:   PetscErrorCode (*vectors)(DS,DSMatType,PetscInt*,PetscReal*);
 31:   PetscErrorCode (*solve[DS_MAX_SOLVE])(DS,PetscScalar*,PetscScalar*);
 32:   PetscErrorCode (*sort)(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
 33:   PetscErrorCode (*sortperm)(DS,PetscInt*,PetscScalar*,PetscScalar*);
 34:   PetscErrorCode (*gettruncatesize)(DS,PetscInt,PetscInt,PetscInt*);
 35:   PetscErrorCode (*truncate)(DS,PetscInt,PetscBool);
 36:   PetscErrorCode (*update)(DS);
 37:   PetscErrorCode (*cond)(DS,PetscReal*);
 38:   PetscErrorCode (*transharm)(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
 39:   PetscErrorCode (*transrks)(DS,PetscScalar);
 40:   PetscErrorCode (*destroy)(DS);
 41:   PetscErrorCode (*matgetsize)(DS,DSMatType,PetscInt*,PetscInt*);
 42:   PetscErrorCode (*hermitian)(DS,DSMatType,PetscBool*);
 43:   PetscErrorCode (*synchronize)(DS,PetscScalar*,PetscScalar*);
 44: };

 46: struct _p_DS {
 47:   PETSCHEADER(struct _DSOps);
 48:   /*------------------------- User parameters --------------------------*/
 49:   DSStateType    state;              /* the current state */
 50:   PetscInt       method;             /* identifies the variant to be used */
 51:   PetscBool      compact;            /* whether the matrices are stored in compact form */
 52:   PetscBool      refined;            /* get refined vectors instead of regular vectors */
 53:   PetscBool      extrarow;           /* assume the matrix dimension is (n+1) x n */
 54:   PetscInt       ld;                 /* leading dimension */
 55:   PetscInt       l;                  /* number of locked (inactive) leading columns */
 56:   PetscInt       n;                  /* current dimension */
 57:   PetscInt       k;                  /* intermediate dimension (e.g. position of arrow) */
 58:   PetscInt       t;                  /* length of decomposition when it was truncated */
 59:   PetscInt       bs;                 /* block size */
 60:   SlepcSC        sc;                 /* sorting criterion */
 61:   DSParallelType pmode;              /* parallel mode (redundant, synchronized, distributed) */

 63:   /*----------------- Status variables and working data ----------------*/
 64:   Mat            omat[DS_NUM_MAT];   /* the matrices (PETSc object) */
 65:   PetscInt       *perm;              /* permutation */
 66:   void           *data;              /* placeholder for solver-specific stuff */
 67:   PetscBool      scset;              /* the sc was provided by the user */
 68:   PetscScalar    *work;
 69:   PetscReal      *rwork;
 70:   PetscBLASInt   *iwork;
 71:   PetscInt       lwork,lrwork,liwork;
 72: };

 74: /*
 75:     Macros to test valid DS arguments
 76: */
 77: #if !defined(PETSC_USE_DEBUG)

 79: #define DSCheckAlloc(h,arg) do {(void)(h);} while (0)
 80: #define DSCheckSolved(h,arg) do {(void)(h);} while (0)
 81: #define DSCheckValidMat(ds,m,arg) do {(void)(ds);} while (0)
 82: #define DSCheckValidMatReal(ds,m,arg) do {(void)(ds);} while (0)

 84: #else

 86: #define DSCheckAlloc(h,arg) \
 87:   do { \
 89:   } while (0)

 91: #define DSCheckSolved(h,arg) \
 92:   do { \
 94:   } while (0)

 96: #define DSCheckValidMat(ds,m,arg) \
 97:   do { \
100:   } while (0)

102: #define DSCheckValidMatReal(ds,m,arg) \
103:   do { \
106:   } while (0)

108: #endif

110: SLEPC_INTERN PetscErrorCode DSAllocateMat_Private(DS,DSMatType);
111: SLEPC_INTERN PetscErrorCode DSAllocateWork_Private(DS,PetscInt,PetscInt,PetscInt);
112: SLEPC_INTERN PetscErrorCode DSSortEigenvalues_Private(DS,PetscScalar*,PetscScalar*,PetscInt*,PetscBool);
113: SLEPC_INTERN PetscErrorCode DSSortEigenvaluesReal_Private(DS,PetscReal*,PetscInt*);
114: SLEPC_INTERN PetscErrorCode DSPermuteColumns_Private(DS,PetscInt,PetscInt,PetscInt,DSMatType,PetscInt*);
115: SLEPC_INTERN PetscErrorCode DSPermuteColumnsTwo_Private(DS,PetscInt,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
116: SLEPC_INTERN PetscErrorCode DSPermuteRows_Private(DS,PetscInt,PetscInt,PetscInt,DSMatType,PetscInt*);
117: SLEPC_INTERN PetscErrorCode DSPermuteBoth_Private(DS,PetscInt,PetscInt,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
118: SLEPC_INTERN PetscErrorCode DSGetTruncateSize_Default(DS,PetscInt,PetscInt,PetscInt*);

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 DSSolve_NHEP_Private(DS,DSMatType,DSMatType,PetscScalar*,PetscScalar*);
129: SLEPC_INTERN PetscErrorCode DSSort_NHEP_Total(DS,DSMatType,DSMatType,PetscScalar*,PetscScalar*);
130: SLEPC_INTERN PetscErrorCode DSSortWithPermutation_NHEP_Private(DS,PetscInt*,DSMatType,DSMatType,PetscScalar*,PetscScalar*);

132: 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);
133: 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);
134: 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);
135: 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);
136: 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*);

138: #endif