Actual source code: dsimpl.h

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2013, 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 PetscLogEvent DS_Solve,DS_Function,DS_Vectors,DS_Other;
 29: PETSC_INTERN const char *DSMatName[];

 31: typedef struct _DSOps *DSOps;

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

 48: struct _p_DS {
 49:   PETSCHEADER(struct _DSOps);
 50:   PetscInt       method;             /* identifies the variant to be used */
 51:   PetscInt       funmethod;          /* to choose among methods for function evaluation */
 52:   PetscBool      compact;            /* whether the matrices are stored in compact form */
 53:   PetscBool      refined;            /* get refined vectors instead of regular vectors */
 54:   PetscBool      extrarow;           /* assume the matrix dimension is (n+1) x n */
 55:   PetscInt       ld;                 /* leading dimension */
 56:   PetscInt       l;                  /* number of locked (inactive) leading columns */
 57:   PetscInt       n;                  /* current dimension */
 58:   PetscInt       m;                  /* current column dimension (for SVD only) */
 59:   PetscInt       k;                  /* intermediate dimension (e.g. position of arrow) */
 60:   PetscInt       t;                  /* length of decomposition when it was truncated */
 61:   DSStateType    state;              /* the current state */
 62:   PetscScalar    *mat[DS_NUM_MAT];   /* the matrices */
 63:   PetscReal      *rmat[DS_NUM_MAT];  /* the matrices (real) */
 64:   PetscInt       *perm;              /* permutation */
 65:   PetscInt       nf;                 /* number of functions in f[] */
 66:   FN             f[DS_NUM_EXTRA];    /* functions provided via DSSetFN() */
 67:   PetscScalar    *work;
 68:   PetscReal      *rwork;
 69:   PetscBLASInt   *iwork;
 70:   PetscInt       lwork,lrwork,liwork;
 71:   /*-------------- User-provided functions and contexts -----------------*/
 72:   PetscErrorCode (*comparison)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 73:   void           *comparisonctx;
 74: };

 76: PETSC_INTERN PetscErrorCode DSAllocateMat_Private(DS,DSMatType);
 77: PETSC_INTERN PetscErrorCode DSAllocateMatReal_Private(DS,DSMatType);
 78: PETSC_INTERN PetscErrorCode DSAllocateWork_Private(DS,PetscInt,PetscInt,PetscInt);
 79: PETSC_INTERN PetscErrorCode DSViewMat_Private(DS,PetscViewer,DSMatType);
 80: PETSC_INTERN PetscErrorCode DSSortEigenvalues_Private(DS,PetscScalar*,PetscScalar*,PetscInt*,PetscBool);
 81: PETSC_INTERN PetscErrorCode DSSortEigenvaluesReal_Private(DS,PetscReal*,PetscInt*);
 82: PETSC_INTERN PetscErrorCode DSPermuteColumns_Private(DS,PetscInt,PetscInt,DSMatType,PetscInt*);
 83: PETSC_INTERN PetscErrorCode DSPermuteRows_Private(DS,PetscInt,PetscInt,DSMatType,PetscInt*);
 84: PETSC_INTERN PetscErrorCode DSPermuteBoth_Private(DS,PetscInt,PetscInt,DSMatType,DSMatType,PetscInt*);
 85: PETSC_INTERN PetscErrorCode DSCopyMatrix_Private(DS,DSMatType,DSMatType);
 86: PETSC_INTERN PetscErrorCode DSSetIdentity(DS,DSMatType);
 87: PETSC_INTERN PetscErrorCode DSComputeMatrix(DS,PetscScalar,PetscBool,DSMatType);
 88: PETSC_INTERN PetscErrorCode DSOrthogonalize(DS,DSMatType,PetscInt,PetscInt*);
 89: PETSC_INTERN PetscErrorCode DSPseudoOrthogonalize(DS,DSMatType,PetscInt,PetscReal*,PetscInt*,PetscReal*);

 91: PETSC_INTERN PetscErrorCode DSGHIEPOrthogEigenv(DS,DSMatType,PetscScalar*,PetscScalar*,PetscBool);
 92: PETSC_INTERN PetscErrorCode DSGHIEPComplexEigs(DS,PetscInt,PetscInt,PetscScalar*,PetscScalar*);
 93: PETSC_INTERN PetscErrorCode DSGHIEPInverseIteration(DS,PetscScalar*,PetscScalar*);
 94: PETSC_INTERN PetscErrorCode DSIntermediate_GHIEP(DS);
 95: PETSC_INTERN PetscErrorCode DSSwitchFormat_GHIEP(DS,PetscBool);
 96: PETSC_INTERN PetscErrorCode DSGHIEPRealBlocks(DS);

 98: PETSC_INTERN PetscErrorCode DSSolve_GHIEP_HZ(DS,PetscScalar*,PetscScalar*);
 99: PETSC_INTERN PetscErrorCode DSSolve_GHIEP_DQDS_II(DS,PetscScalar*,PetscScalar*);

101: #endif