Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1908 jroman 1
/*                      
2
 
3
   Private header for QEPLINEAR.
4
 
5
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 7
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
1908 jroman 8
 
9
   This file is part of SLEPc.
10
 
11
   SLEPc is free software: you can redistribute it and/or modify it under  the
12
   terms of version 3 of the GNU Lesser General Public License as published by
13
   the Free Software Foundation.
14
 
15
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
16
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
17
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
18
   more details.
19
 
20
   You  should have received a copy of the GNU Lesser General  Public  License
21
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23
*/
24
 
25
typedef struct {
2216 jroman 26
  PetscBool  explicitmatrix;
1909 jroman 27
  PetscInt   cform;           /* companion form */
1919 jroman 28
  PetscReal  sfactor;         /* scaling factor */
1908 jroman 29
  Mat        A,B;             /* matrices of generalized eigenproblem */
30
  EPS        eps;             /* linear eigensolver for Az=lBz */
31
  Mat        M,C,K;           /* copy of QEP coefficient matrices */
32
  Vec        x1,x2,y1,y2;     /* work vectors */
2216 jroman 33
  PetscBool  setfromoptionscalled;
1908 jroman 34
} QEP_LINEAR;
35
 
1910 jroman 36
/* N1 */
2324 jroman 37
extern PetscErrorCode MatMult_Linear_N1A(Mat,Vec,Vec);
38
extern PetscErrorCode MatMult_Linear_N1B(Mat,Vec,Vec);
39
extern PetscErrorCode MatGetDiagonal_Linear_N1A(Mat,Vec);
40
extern PetscErrorCode MatGetDiagonal_Linear_N1B(Mat,Vec);
41
extern PetscErrorCode MatCreateExplicit_Linear_N1A(MPI_Comm,QEP_LINEAR*,Mat*);
42
extern PetscErrorCode MatCreateExplicit_Linear_N1B(MPI_Comm,QEP_LINEAR*,Mat*);
1908 jroman 43
 
1910 jroman 44
/* N2 */
2324 jroman 45
extern PetscErrorCode MatMult_Linear_N2A(Mat,Vec,Vec);
46
extern PetscErrorCode MatMult_Linear_N2B(Mat,Vec,Vec);
47
extern PetscErrorCode MatGetDiagonal_Linear_N2A(Mat,Vec);
48
extern PetscErrorCode MatGetDiagonal_Linear_N2B(Mat,Vec);
49
extern PetscErrorCode MatCreateExplicit_Linear_N2A(MPI_Comm,QEP_LINEAR*,Mat*);
50
extern PetscErrorCode MatCreateExplicit_Linear_N2B(MPI_Comm,QEP_LINEAR*,Mat*);
1910 jroman 51
 
1911 jroman 52
/* S1 */
2324 jroman 53
extern PetscErrorCode MatMult_Linear_S1A(Mat,Vec,Vec);
54
extern PetscErrorCode MatMult_Linear_S1B(Mat,Vec,Vec);
55
extern PetscErrorCode MatGetDiagonal_Linear_S1A(Mat,Vec);
56
extern PetscErrorCode MatGetDiagonal_Linear_S1B(Mat,Vec);
57
extern PetscErrorCode MatCreateExplicit_Linear_S1A(MPI_Comm,QEP_LINEAR*,Mat*);
58
extern PetscErrorCode MatCreateExplicit_Linear_S1B(MPI_Comm,QEP_LINEAR*,Mat*);
1911 jroman 59
 
1912 jroman 60
/* S2 */
2324 jroman 61
extern PetscErrorCode MatMult_Linear_S2A(Mat,Vec,Vec);
62
extern PetscErrorCode MatMult_Linear_S2B(Mat,Vec,Vec);
63
extern PetscErrorCode MatGetDiagonal_Linear_S2A(Mat,Vec);
64
extern PetscErrorCode MatGetDiagonal_Linear_S2B(Mat,Vec);
65
extern PetscErrorCode MatCreateExplicit_Linear_S2A(MPI_Comm,QEP_LINEAR*,Mat*);
66
extern PetscErrorCode MatCreateExplicit_Linear_S2B(MPI_Comm,QEP_LINEAR*,Mat*);
1912 jroman 67
 
1913 jroman 68
/* H1 */
2324 jroman 69
extern PetscErrorCode MatMult_Linear_H1A(Mat,Vec,Vec);
70
extern PetscErrorCode MatMult_Linear_H1B(Mat,Vec,Vec);
71
extern PetscErrorCode MatGetDiagonal_Linear_H1A(Mat,Vec);
72
extern PetscErrorCode MatGetDiagonal_Linear_H1B(Mat,Vec);
73
extern PetscErrorCode MatCreateExplicit_Linear_H1A(MPI_Comm,QEP_LINEAR*,Mat*);
74
extern PetscErrorCode MatCreateExplicit_Linear_H1B(MPI_Comm,QEP_LINEAR*,Mat*);
1913 jroman 75
 
1914 jroman 76
/* H2 */
2324 jroman 77
extern PetscErrorCode MatMult_Linear_H2A(Mat,Vec,Vec);
78
extern PetscErrorCode MatMult_Linear_H2B(Mat,Vec,Vec);
79
extern PetscErrorCode MatGetDiagonal_Linear_H2A(Mat,Vec);
80
extern PetscErrorCode MatGetDiagonal_Linear_H2B(Mat,Vec);
81
extern PetscErrorCode MatCreateExplicit_Linear_H2A(MPI_Comm,QEP_LINEAR*,Mat*);
82
extern PetscErrorCode MatCreateExplicit_Linear_H2B(MPI_Comm,QEP_LINEAR*,Mat*);
1914 jroman 83