Actual source code: agmresimpl.h

  1: /*
  2:   Private data structure used for the KSP AGMRES.
  3:   It extends the definition of KSP_GMRES and KSP_DGMRES data structures.
  4:   If you modify something there (located in gmresimpl.h and in dgmresimpl.h), you should  modify it here as well.
  5:   In this KSP, KSPSIZE denotes the size of the basis (possibly augmented with Schur vectors) and MAXKSPSIZE denotes the maximum size of the augmented basis
  6: */
  7: #pragma once

  9: #include <../src/ksp/ksp/impls/gmres/dgmres/dgmresimpl.h>
 10: typedef struct {
 11:   KSPGMRESHEADER
 12:   KSPDGMRESHEADER

 14:   /* Data specific to AGMRES */
 15:   PetscReal     bgv;                        /* large multiple of the remaining allowed number of steps -- used for the adaptive strategy */
 16:   PetscBool     ritz;                       /* Compute the Harmonic Ritz vectors instead of the Ritz vectors */
 17:   PetscBool     DeflPrecond;                /* Apply deflation by building adaptively a preconditioner, otherwise augment the basis */
 18:   PetscScalar  *Qloc;                       /* Orthogonal reflectors from the QR of the basis */
 19:   PetscScalar  *Rloc;                       /* triangular matrix obtained from the QR of the basis */
 20:   PetscScalar  *Rshift, *Ishift;            /* Real and Imaginary parts of the shifts in the Newton basis */
 21:   PetscScalar  *Scale;                      /* Norm of the vectors in the Newton basis */
 22:   PetscBool     HasShifts;                  /* Estimation of shifts exists */
 23:   PetscMPIInt   rank, size;                 /* Rank and size of the current process; to be used in RODDEC*/
 24:   PetscMPIInt   First, Last, Ileft, Iright; /* Create a ring of processors for RODDEC */
 25:   PetscScalar  *MatEigL, *MatEigR;          /* matrices for the eigenvalue problem */
 26:   PetscScalar  *sgn;                        /* Sign of the rotation in the QR factorization of the basis */
 27:   PetscScalar  *tloc;                       /* */
 28:   Vec          *TmpU;                       /* Temporary vectors */
 29:   PetscScalar  *beta;                       /* needed for the eigenvalues */
 30:   PetscBLASInt *select;                     /* array used to select the Schur vectors to order */
 31:   PetscScalar  *temp, *wbufptr;
 32:   PetscScalar  *tau; /* Scalar factors of the elementary reflectors in xgeqrf */
 33:   PetscMPIInt   tag;
 34: } KSP_AGMRES;

 36: PETSC_EXTERN PetscLogEvent KSP_AGMRESComputeDeflationData;
 37: PETSC_EXTERN PetscLogEvent KSP_AGMRESBuildBasis;
 38: PETSC_EXTERN PetscLogEvent KSP_AGMRESComputeShifts;
 39: PETSC_EXTERN PetscLogEvent KSP_AGMRESRoddec;

 41: /* vector names */
 42: #define VEC_TMP       agmres->vecs[0]
 43: #define VEC_TMP_MATOP agmres->vecs[1]
 44: #define VEC_V(i)      agmres->vecs[VEC_OFFSET + i]

 46: #define MAXKSPSIZE ((agmres->DeflPrecond) ? (agmres->max_k) : (agmres->max_k + agmres->max_neig))
 47: #define KSPSIZE    ((agmres->DeflPrecond) ? (agmres->max_k) : (agmres->max_k + agmres->r))
 48: #define H(a, b)    (agmres->hh_origin + (b) * (MAXKSPSIZE + 2) + (a))
 49: #define HS(a, b)   (agmres->hes_origin + (b) * (MAXKSPSIZE + 1) + (a))
 50: #define RLOC(a, b) (agmres->Rloc + (b) * (MAXKSPSIZE + 1) + (a))

 52: PetscErrorCode KSPAGMRESRoddec(KSP, PetscInt);
 53: PetscErrorCode KSPAGMRESRodvec(KSP, PetscInt, PetscScalar *, Vec);
 54: PetscErrorCode KSPAGMRESLejaOrdering(PetscScalar *, PetscScalar *, PetscScalar *, PetscScalar *, PetscInt);
 55: PetscErrorCode KSPAGMRESRoddecInitNeighboor(KSP);
 56: PetscErrorCode KSPAGMRESComputeDeflationData(KSP);

 58: #define AGMRES_DEFAULT_MAXK     30
 59: #define AGMRES_DELTA_DIRECTIONS 10