Actual source code: dgmresimpl.h

petsc-3.7.5 2017-01-01
Report Typos and Errors
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: #define KSPDGMRESHEADER \ 10: /* Data specific to DGMRES */ \ 11: Vec *U; /* Vectors that form the basis of the invariant subspace */ \ 12: PetscScalar *T; /* T=U^T*M^{-1}*A*U */ \ 13: PetscScalar *TF; /* The factors L and U from T = P*L*U */ \ 14: PetscBLASInt *InvP; /* Permutation Vector from the LU factorization of T */ \ 15: PetscInt neig; /* number of eigenvalues to extract at each restart */ \ 16: PetscInt r; /* current number of deflated eigenvalues */ \ 17: PetscInt max_neig; /* Maximum number of eigenvalues to deflate */ \ 18: PetscReal lambdaN; /* modulus of the largest eigenvalue of A */ \ 19: PetscReal smv; /* smaller multiple of the remaining allowed number of steps -- used for the adaptive strategy */ \ 20: PetscBool force; /* Force the use of the deflation at the restart */ \ 21: PetscInt matvecs; /* Total number of matrix-vectors */ \ 22: PetscInt GreatestEig; /* Extract the greatest eigenvalues instead */ \ 23: PetscReal *wr, *wi, *modul; /* Real and complex part and modulus of eigenvalues */ \ 24: PetscScalar *Q, *Z; /* Left and right schur/eigenvectors from the QZ algorithm */ \ 25: PetscInt *perm; /* temporary permutation vector */ \ 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: PetscScalar *work; /* work space for LAPACK functions */ \ 40: PetscBLASInt *iwork; /* work space for LAPACK functions */ \ 41: PetscReal *orth; /* Coefficients for the orthogonalization */ \ 42: PetscBool HasSchur; /* Indicate if the Schur form had already been computed in this cycle */ \ 43: PetscBool improve; /* 0 = do not improve the eigenvalues; This is an experimental option */ 45: typedef struct { 46: KSPGMRESHEADER 47: KSPDGMRESHEADER 48: } KSP_DGMRES; 50: PETSC_EXTERN PetscLogEvent KSP_DGMRESComputeDeflationData, KSP_DGMRESApplyDeflation; 51: #define HH(a,b) (dgmres->hh_origin + (b)*(dgmres->max_k+2)+(a)) 52: #define HES(a,b) (dgmres->hes_origin + (b)*(dgmres->max_k+1)+(a)) 53: #define CC(a) (dgmres->cc_origin + (a)) 54: #define SS(a) (dgmres->ss_origin + (a)) 55: #define GRS(a) (dgmres->rs_origin + (a)) 57: /* vector names */ 58: #define VEC_OFFSET 2 59: #define VEC_TEMP dgmres->vecs[0] 60: #define VEC_TEMP_MATOP dgmres->vecs[1] 61: #define VEC_VV(i) dgmres->vecs[VEC_OFFSET+i] 63: #define EIG_OFFSET 1 64: #define DGMRES_DEFAULT_EIG 1 65: #define DGMRES_DEFAULT_MAXEIG 10 67: #define UU dgmres->U 68: #define TT dgmres->T 69: #define TTF dgmres->TF 70: #define XX dgmres->X 71: #define INVP dgmres->InvP 72: #define MU dgmres->mu 73: #define MX dgmres->mx 74: #define UMX dgmres->umx 75: #define XMU dgmres->xmu 76: #define XMX dgmres->xmx 77: #define X1 dgmres->x1 78: #define X2 dgmres->x2 79: #define SR dgmres->Sr 80: #define SR2 dgmres->Sr2 81: #define AUAU dgmres->auau 82: #define AUU dgmres->auu 83: #define MAX_K dgmres->max_k 84: #define MAX_NEIG dgmres->max_neig 85: #define WORK dgmres->work 86: #define IWORK dgmres->iwork 87: #define ORTH dgmres->orth 88: #define SMV 1 89: #endif