Actual source code: dgmresimpl.h

  1: #pragma once

  3: #define KSPGMRES_NO_MACROS
  4: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
  5: #include <petscblaslapack.h>

  7: #define KSPDGMRESHEADER \
  8:   /* Data specific to DGMRES */ \
  9:   Vec          *U;               /* Vectors that form the basis of the invariant subspace */ \
 10:   PetscScalar  *T;               /* T=U^T*M^{-1}*A*U */ \
 11:   PetscScalar  *TF;              /* The factors L and U from T = P*L*U */ \
 12:   PetscBLASInt *InvP;            /* Permutation Vector from the LU factorization of T */ \
 13:   PetscInt      neig;            /* number of eigenvalues to extract at each restart */ \
 14:   PetscInt      r;               /* current number of deflated eigenvalues */ \
 15:   PetscInt      max_neig;        /* Maximum number of eigenvalues to deflate */ \
 16:   PetscReal     lambdaN;         /* modulus of the largest eigenvalue of A */ \
 17:   PetscReal     smv;             /* smaller multiple of the remaining allowed number of steps -- used for the adaptive strategy */ \
 18:   PetscBool     force;           /* Force the use of the deflation at the restart */ \
 19:   PetscInt      matvecs;         /* Total number of matrix-vectors */ \
 20:   PetscInt      GreatestEig;     /* Extract the greatest eigenvalues instead */ \
 21:   PetscReal    *wr, *wi, *modul; /* Real and complex part and modulus of eigenvalues */ \
 22:   PetscScalar  *Q, *Z;           /* Left and right schur/eigenvectors from the QZ algorithm */ \
 23:   PetscInt     *perm;            /* temporary permutation vector */ \
 24:   /* Work spaces */ \
 25:   Vec          *mu;       /* Save the product M^{-1}AU */ \
 26:   PetscScalar  *Sr;       /* Schur vectors to extract */ \
 27:   Vec          *X;        /* Schurs Vectors Sr projected to the entire space */ \
 28:   Vec          *mx;       /* store the product M^{-1}*A*X */ \
 29:   PetscScalar  *umx;      /* store the product U^T*M^{-1}*A*X */ \
 30:   PetscScalar  *xmu;      /* store the product X^T*M^{-1}*A*U */ \
 31:   PetscScalar  *xmx;      /* store the product X^T*M^{-1}*A*X */ \
 32:   PetscScalar  *x1;       /* store the product U^T*x */ \
 33:   PetscScalar  *x2;       /* store the product U^T*x */ \
 34:   PetscScalar  *Sr2;      /* Schur vectors at the improvement step */ \
 35:   PetscScalar  *auau;     /* product of (M*A*U)^T*M*A*U */ \
 36:   PetscScalar  *auu;      /* product of (M*A*U)^T*U */ \
 37:   PetscScalar  *work;     /* work space for LAPACK functions */ \
 38:   PetscBLASInt *iwork;    /* work space for LAPACK functions */ \
 39:   PetscReal    *orth;     /* Coefficients for the orthogonalization */ \
 40:   PetscBool     HasSchur; /* Indicate if the Schur form had already been computed in this cycle */ \
 41:   PetscBool     improve;  /* 0 = do not improve the eigenvalues; This is an experimental option */

 43: typedef struct {
 44:   KSPGMRESHEADER
 45:   KSPDGMRESHEADER
 46: } KSP_DGMRES;

 48: PETSC_INTERN PetscErrorCode KSPDGMRESComputeDeflationData(KSP, PetscInt *);

 50: PETSC_EXTERN PetscLogEvent KSP_DGMRESComputeDeflationData;
 51: PETSC_EXTERN PetscLogEvent KSP_DGMRESApplyDeflation;

 53: #define HH(a, b)  (dgmres->hh_origin + (b) * (dgmres->max_k + 2) + (a))
 54: #define HES(a, b) (dgmres->hes_origin + (b) * (dgmres->max_k + 1) + (a))
 55: #define CC(a)     (dgmres->cc_origin + (a))
 56: #define SS(a)     (dgmres->ss_origin + (a))
 57: #define GRS(a)    (dgmres->rs_origin + (a))

 59: /* vector names */
 60: #define VEC_OFFSET     2
 61: #define VEC_TEMP       dgmres->vecs[0]
 62: #define VEC_TEMP_MATOP dgmres->vecs[1]
 63: #define VEC_VV(i)      dgmres->vecs[VEC_OFFSET + i]

 65: #define EIG_OFFSET            1
 66: #define DGMRES_DEFAULT_EIG    1
 67: #define DGMRES_DEFAULT_MAXEIG 10

 69: #define UU       dgmres->U
 70: #define TT       dgmres->T
 71: #define TTF      dgmres->TF
 72: #define XX       dgmres->X
 73: #define INVP     dgmres->InvP
 74: #define MU       dgmres->mu
 75: #define MX       dgmres->mx
 76: #define UMX      dgmres->umx
 77: #define XMU      dgmres->xmu
 78: #define XMX      dgmres->xmx
 79: #define X1       dgmres->x1
 80: #define X2       dgmres->x2
 81: #define SR       dgmres->Sr
 82: #define SR2      dgmres->Sr2
 83: #define AUAU     dgmres->auau
 84: #define AUU      dgmres->auu
 85: #define MAX_K    dgmres->max_k
 86: #define MAX_NEIG dgmres->max_neig
 87: #define WORK     dgmres->work
 88: #define IWORK    dgmres->iwork
 89: #define ORTH     dgmres->orth
 90: #define SMV      1

 92: #define GMRES_DELTA_DIRECTIONS 10
 93: #define GMRES_DEFAULT_MAXK     30