Actual source code: dgmresimpl.h
petsc-3.3-p5 2012-12-01
4: #include <petsc-private/kspimpl.h> /*I "petscksp.h" I*/
5: #include <petscblaslapack.h>
6: #define KSPGMRES_NO_MACROS
7: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
9: typedef struct {
10: KSPGMRESHEADER
12: /* Data specific to DGMRES */
13: Vec *U; /* Vectors that form the basis of the invariant subspace */
14: PetscScalar *T; /* T=U^T*M^{-1}*A*U */
15: PetscScalar *TF; /* The factors L and U from T = P*L*U */
16: PetscBLASInt *InvP; /* Permutation Vector from the LU factorization of T */
17: PetscInt neig; /* number of eigenvalues to extract at each restart */
18: PetscInt r; /* current number of deflated eigenvalues */
19: PetscInt max_neig; /* Maximum number of eigenvalues to deflate */
20: PetscReal lambdaN; /* modulus of the largest eigenvalue of A */
21: PetscReal smv; /* smaller multiple of the remaining allowed number of steps -- used for the adaptive strategy */
22: PetscInt force; /* Force the use of the deflation at the restart */
23: PetscInt matvecs; /* Total number of matrix-vectors */
24: PetscInt GreatestEig; /* Extract the greatest eigenvalues instead */
25:
26: /* Work spaces */
27: Vec *mu; /* Save the product M^{-1}AU */
28: PetscScalar *Sr; /* Schur vectors to extract */
29: Vec *X; /* Schurs Vectors Sr projected to the entire space */
30: Vec *mx; /* store the product M^{-1}*A*X */
31: PetscScalar *umx; /* store the product U^T*M^{-1}*A*X */
32: PetscScalar *xmu; /* store the product X^T*M^{-1}*A*U */
33: PetscScalar *xmx; /* store the product X^T*M^{-1}*A*X */
34: PetscScalar *x1; /* store the product U^T*x */
35: PetscScalar *x2; /* store the product U^T*x */
36: PetscScalar *Sr2; /* Schur vectors at the improvement step */
37: PetscScalar *auau; /* product of (M*A*U)^T*M*A*U */
38: PetscScalar *auu; /* product of (M*A*U)^T*U */
39:
40: PetscScalar *work; /* work space for LAPACK functions */
41: PetscInt *iwork; /* work space for LAPACK functions */
42: PetscReal *orth; /* Coefficients for the orthogonalization */
43:
44: PetscInt improve; /* 0 = do not improve the eigenvalues; This is an experimental option */
45:
46: } KSP_DGMRES;
48: PetscLogEvent KSP_DGMRESComputeDeflationData, KSP_DGMRESApplyDeflation;
49: #define HH(a,b) (dgmres->hh_origin + (b)*(dgmres->max_k+2)+(a))
50: #define HES(a,b) (dgmres->hes_origin + (b)*(dgmres->max_k+1)+(a))
51: #define CC(a) (dgmres->cc_origin + (a))
52: #define SS(a) (dgmres->ss_origin + (a))
53: #define GRS(a) (dgmres->rs_origin + (a))
55: /* vector names */
56: #define VEC_OFFSET 2
57: #define VEC_TEMP dgmres->vecs[0]
58: #define VEC_TEMP_MATOP dgmres->vecs[1]
59: #define VEC_VV(i) dgmres->vecs[VEC_OFFSET+i]
61: #define EIG_OFFSET 2
62: #define DGMRES_DEFAULT_EIG 1
63: #define DGMRES_DEFAULT_MAXEIG 100
65: #define UU dgmres->U
66: #define TT dgmres->T
67: #define TTF dgmres->TF
68: #define XX dgmres->X
69: #define INVP dgmres->InvP
70: #define MU dgmres->mu
71: #define MX dgmres->mx
72: #define UMX dgmres->umx
73: #define XMU dgmres->xmu
74: #define XMX dgmres->xmx
75: #define X1 dgmres->x1
76: #define X2 dgmres->x2
77: #define SR dgmres->Sr
78: #define SR2 dgmres->Sr2
79: #define AUAU dgmres->auau
80: #define AUU dgmres->auu
81: #define MAX_K dgmres->max_k
82: #define MAX_NEIG dgmres->max_neig
83: #define WORK dgmres->work
84: #define IWORK dgmres->iwork
85: #define ORTH dgmres->orth
86: #define SMV 1
87: #endif