1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) , 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(_STIMPL)
23: #define _STIMPL 25: #include <slepcst.h>
26: #include <slepc/private/slepcimpl.h>
28: PETSC_EXTERN PetscBool STRegisterAllCalled;
29: PETSC_EXTERN PetscErrorCode STRegisterAll(void);
30: PETSC_EXTERN PetscLogEvent ST_SetUp,ST_Apply,ST_ApplyTranspose,ST_MatSetUp,ST_MatMult,ST_MatMultTranspose,ST_MatSolve,ST_MatSolveTranspose;
32: typedef struct _STOps *STOps;
34: struct _STOps {
35: PetscErrorCode (*setup)(ST);
36: PetscErrorCode (*apply)(ST,Vec,Vec);
37: PetscErrorCode (*getbilinearform)(ST,Mat*);
38: PetscErrorCode (*applytrans)(ST,Vec,Vec);
39: PetscErrorCode (*setshift)(ST,PetscScalar);
40: PetscErrorCode (*setfromoptions)(PetscOptionItems*,ST);
41: PetscErrorCode (*postsolve)(ST);
42: PetscErrorCode (*backtransform)(ST,PetscInt,PetscScalar*,PetscScalar*);
43: PetscErrorCode (*destroy)(ST);
44: PetscErrorCode (*reset)(ST);
45: PetscErrorCode (*view)(ST,PetscViewer);
46: PetscErrorCode (*checknullspace)(ST,BV);
47: };
49: /*
50: 'Updated' state means STSetUp must be called because matrices have been
51: modified, but the pattern is the same (hence reuse symbolic factorization)
52: */
53: typedef enum { ST_STATE_INITIAL,
54: ST_STATE_SETUP,
55: ST_STATE_UPDATED } STStateType;
57: struct _p_ST {
58: PETSCHEADER(struct _STOps);
59: /*------------------------- User parameters --------------------------*/
60: Mat *A; /* matrices that define the eigensystem */
61: PetscObjectState *Astate; /* state (to identify the original matrices) */
62: Mat *T; /* matrices resulting from transformation */
63: Mat P; /* matrix from which preconditioner is built */
64: PetscInt nmat; /* number of matrices */
65: PetscScalar sigma; /* value of the shift */
66: PetscBool sigma_set; /* whether the user provided the shift or not */
67: PetscScalar defsigma; /* default value of the shift */
68: STMatMode shift_matrix;
69: MatStructure str; /* whether matrices have the same pattern or not */
70: PetscBool transform; /* whether transformed matrices are computed */
72: /*------------------------- Misc data --------------------------*/
73: KSP ksp;
74: Vec w; /* work vector used in apply operation */
75: Vec D; /* diagonal matrix for balancing */
76: Vec wb; /* balancing requires an extra work vector */
77: void *data;
78: STStateType state; /* initial -> setup -> with updated matrices */
79: };
83: /*
84: ST_AllocateWorkVec - Allocate work vector for the STApply operation.
85: */
86: PETSC_STATIC_INLINE PetscErrorCode ST_AllocateWorkVec(ST st) 87: {
91: if (!st->w) {
92: MatCreateVecs(st->A[0],&st->w,NULL);
93: PetscLogObjectParent((PetscObject)st,(PetscObject)st->w);
94: }
95: return(0);
96: }
98: /*
99: Macros to test valid ST arguments
100: */
101: #if !defined(PETSC_USE_DEBUG)
103: #define STCheckMatrices(h,arg) do {} while (0)105: #else
107: #define STCheckMatrices(h,arg) \108: do { \109: if (!h->A) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"ST matrices have not been set: Parameter #%d",arg); \110: } while (0)112: #endif
114: PETSC_INTERN PetscErrorCode STGetBilinearForm_Default(ST,Mat*);
115: PETSC_INTERN PetscErrorCode STCheckNullSpace_Default(ST,BV);
116: PETSC_INTERN PetscErrorCode STMatShellCreate(ST,PetscScalar,PetscInt,PetscInt*,PetscScalar*,Mat*);
117: PETSC_INTERN PetscErrorCode STMatShellShift(Mat,PetscScalar);
118: PETSC_INTERN PetscErrorCode STMatSetHermitian(ST,Mat);
119: PETSC_INTERN PetscErrorCode STCheckFactorPackage(ST);
120: PETSC_INTERN PetscErrorCode STMatMAXPY_Private(ST,PetscScalar,PetscScalar,PetscInt,PetscScalar*,PetscBool,Mat*);
121: PETSC_INTERN PetscErrorCode STCoeffs_Monomial(ST,PetscScalar*);
123: #endif