Actual source code: dsimpl.h

slepc-3.20.2 2024-03-15
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: #pragma once

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

 16: /* SUBMANSEC = DS */

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

 23: typedef struct _DSOps *DSOps;

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

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

 62:   /*----------------- Status variables and working data ----------------*/
 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 {(void)(h);} while (0)
 79: #define DSCheckSolved(h,arg) do {(void)(h);} while (0)
 80: #define DSCheckValidMat(ds,m,arg) do {(void)(ds);} while (0)
 81: #define DSCheckValidMatReal(ds,m,arg) do {(void)(ds);} while (0)

 83: #else

 85: #define DSCheckAlloc(h,arg) \
 86:   do { \
 87:     PetscCheck((h)->ld,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:     PetscCheck((h)->state>=DS_STATE_CONDENSED,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:     PetscCheck((m)<DS_NUM_MAT,PetscObjectComm((PetscObject)(ds)),PETSC_ERR_ARG_WRONG,"Invalid matrix: Parameter #%d",arg); \
 98:     PetscCheck((ds)->omat[m],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:     PetscCheck((m)==DS_MAT_T || (m)==DS_MAT_D,PetscObjectComm((PetscObject)(ds)),PETSC_ERR_ARG_WRONG,"Invalid matrix, can only be used for T and D: Parameter #%d",arg); \
104:     PetscCheck((ds)->omat[m],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 DSAllocateMat_Private(DS,DSMatType);
110: SLEPC_INTERN PetscErrorCode DSAllocateWork_Private(DS,PetscInt,PetscInt,PetscInt);
111: SLEPC_INTERN PetscErrorCode DSSortEigenvalues_Private(DS,PetscScalar*,PetscScalar*,PetscInt*,PetscBool);
112: SLEPC_INTERN PetscErrorCode DSSortEigenvaluesReal_Private(DS,PetscReal*,PetscInt*);
113: SLEPC_INTERN PetscErrorCode DSPermuteColumns_Private(DS,PetscInt,PetscInt,PetscInt,DSMatType,PetscInt*);
114: SLEPC_INTERN PetscErrorCode DSPermuteColumnsTwo_Private(DS,PetscInt,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
115: SLEPC_INTERN PetscErrorCode DSPermuteRows_Private(DS,PetscInt,PetscInt,PetscInt,DSMatType,PetscInt*);
116: SLEPC_INTERN PetscErrorCode DSPermuteBoth_Private(DS,PetscInt,PetscInt,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
117: SLEPC_INTERN PetscErrorCode DSGetTruncateSize_Default(DS,PetscInt,PetscInt,PetscInt*);

119: SLEPC_INTERN PetscErrorCode DSGHIEPOrthogEigenv(DS,DSMatType,PetscScalar*,PetscScalar*,PetscBool);
120: SLEPC_INTERN PetscErrorCode DSPseudoOrthog_HR(PetscInt*,PetscScalar*,PetscInt,PetscReal*,PetscScalar*,PetscInt,PetscBLASInt*,PetscBLASInt*,PetscBool*,PetscScalar*);
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*);
127: SLEPC_INTERN PetscErrorCode DSArrowTridiag(PetscBLASInt,PetscReal*,PetscReal*,PetscScalar*,PetscBLASInt);

129: SLEPC_INTERN PetscErrorCode DSSolve_NHEP_Private(DS,DSMatType,DSMatType,PetscScalar*,PetscScalar*);
130: SLEPC_INTERN PetscErrorCode DSSort_NHEP_Total(DS,DSMatType,DSMatType,PetscScalar*,PetscScalar*);
131: SLEPC_INTERN PetscErrorCode DSSortWithPermutation_NHEP_Private(DS,PetscInt*,DSMatType,DSMatType,PetscScalar*,PetscScalar*);

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