Actual source code: dsimpl.h

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: #if !defined(_DSIMPL)
 23: #define _DSIMPL

 25: #include <slepcds.h>
 26: #include <slepc/private/slepcimpl.h>

 28: PETSC_EXTERN PetscBool DSRegisterAllCalled;
 29: PETSC_EXTERN PetscErrorCode DSRegisterAll(void);
 30: PETSC_EXTERN PetscLogEvent DS_Solve,DS_Vectors,DS_Other;
 31: PETSC_INTERN const char *DSMatName[];

 33: typedef struct _DSOps *DSOps;

 35: struct _DSOps {
 36:   PetscErrorCode (*allocate)(DS,PetscInt);
 37:   PetscErrorCode (*view)(DS,PetscViewer);
 38:   PetscErrorCode (*vectors)(DS,DSMatType,PetscInt*,PetscReal*);
 39:   PetscErrorCode (*solve[DS_MAX_SOLVE])(DS,PetscScalar*,PetscScalar*);
 40:   PetscErrorCode (*sort)(DS,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*);
 41:   PetscErrorCode (*truncate)(DS,PetscInt);
 42:   PetscErrorCode (*update)(DS);
 43:   PetscErrorCode (*cond)(DS,PetscReal*);
 44:   PetscErrorCode (*transharm)(DS,PetscScalar,PetscReal,PetscBool,PetscScalar*,PetscReal*);
 45:   PetscErrorCode (*transrks)(DS,PetscScalar);
 46:   PetscErrorCode (*normalize)(DS,DSMatType,PetscInt);
 47:   PetscErrorCode (*destroy)(DS);
 48: };

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

 67:   /*----------------- Status variables and working data ----------------*/
 68:   PetscScalar    *mat[DS_NUM_MAT];   /* the matrices */
 69:   PetscReal      *rmat[DS_NUM_MAT];  /* the matrices (real) */
 70:   Mat            omat[DS_NUM_MAT];   /* the matrices (PETSc object) */
 71:   PetscInt       *perm;              /* permutation */
 72:   void           *data;              /* placeholder for solver-specific stuff */
 73:   PetscScalar    *work;
 74:   PetscReal      *rwork;
 75:   PetscBLASInt   *iwork;
 76:   PetscInt       lwork,lrwork,liwork;
 77: };

 79: /*
 80:     Macros to test valid DS arguments
 81: */
 82: #if !defined(PETSC_USE_DEBUG)

 84: #define DSCheckAlloc(h,arg) do {} while (0)
 85: #define DSCheckSolved(h,arg) do {} while (0)

 87: #else

 89: #define DSCheckAlloc(h,arg) \
 90:   do { \
 91:     if (!h->ld) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call DSAllocate() first: Parameter #%d",arg); \
 92:   } while (0)

 94: #define DSCheckSolved(h,arg) \
 95:   do { \
 96:     if (h->state<DS_STATE_CONDENSED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call DSSolve() first: Parameter #%d",arg); \
 97:   } while (0)

 99: #endif

101: PETSC_INTERN PetscErrorCode DSAllocateMat_Private(DS,DSMatType);
102: PETSC_INTERN PetscErrorCode DSAllocateMatReal_Private(DS,DSMatType);
103: PETSC_INTERN PetscErrorCode DSAllocateWork_Private(DS,PetscInt,PetscInt,PetscInt);
104: PETSC_INTERN PetscErrorCode DSSortEigenvalues_Private(DS,PetscScalar*,PetscScalar*,PetscInt*,PetscBool);
105: PETSC_INTERN PetscErrorCode DSSortEigenvaluesReal_Private(DS,PetscReal*,PetscInt*);
106: PETSC_INTERN PetscErrorCode DSPermuteColumns_Private(DS,PetscInt,PetscInt,DSMatType,PetscInt*);
107: PETSC_INTERN PetscErrorCode DSPermuteRows_Private(DS,PetscInt,PetscInt,DSMatType,PetscInt*);
108: PETSC_INTERN PetscErrorCode DSPermuteBoth_Private(DS,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
109: PETSC_INTERN PetscErrorCode DSCopyMatrix_Private(DS,DSMatType,DSMatType);

111: PETSC_INTERN PetscErrorCode DSGHIEPOrthogEigenv(DS,DSMatType,PetscScalar*,PetscScalar*,PetscBool);
112: PETSC_INTERN PetscErrorCode DSGHIEPComplexEigs(DS,PetscInt,PetscInt,PetscScalar*,PetscScalar*);
113: PETSC_INTERN PetscErrorCode DSGHIEPInverseIteration(DS,PetscScalar*,PetscScalar*);
114: PETSC_INTERN PetscErrorCode DSIntermediate_GHIEP(DS);
115: PETSC_INTERN PetscErrorCode DSSwitchFormat_GHIEP(DS,PetscBool);
116: PETSC_INTERN PetscErrorCode DSGHIEPRealBlocks(DS);

118: PETSC_INTERN PetscErrorCode DSSolve_GHIEP_HZ(DS,PetscScalar*,PetscScalar*);
119: PETSC_INTERN PetscErrorCode DSSolve_GHIEP_DQDS_II(DS,PetscScalar*,PetscScalar*);

121: 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);
122: 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);
123: 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);
124: 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);
125: 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*);

127: #endif