Actual source code: dgmresimpl.h

  1: /*
  2:    Private data structure used by the GMRES method. This data structure
  3:   must be identical to the beginning of the KSP_FGMRES data structure
  4:   so if you CHANGE anything here you must also change it there.
  5: */

  9: #include <private/kspimpl.h>        /*I "petscksp.h" I*/
 10: #include <petscblaslapack.h>

 12: typedef struct {
 13:   /* Hessenberg matrix and orthogonalization information.  Hes holds
 14:        the original (unmodified) hessenberg matrix which may be used
 15:        to estimate the Singular Values of the matrix */
 16:   PetscScalar *hh_origin,*hes_origin,*cc_origin,*ss_origin,*rs_origin;

 18:   PetscScalar *orthogwork; /* holds dot products computed in orthogonalization */

 20:   /* Work space for computing eigenvalues/singular values */
 21:   PetscReal   *Dsvd;
 22:   PetscScalar *Rsvd;
 23: 
 24:   /* parameters */
 25:   PetscReal   haptol;            /* tolerance used for the "HAPPY BREAK DOWN"  */
 26:   PetscInt    max_k;             /* maximum size of the approximation space  
 27:   before restarting */

 29:   PetscErrorCode (*orthog)(KSP,PetscInt); /* Functions to use (special to gmres) */
 30:   KSPGMRESCGSRefinementType cgstype;
 31: 
 32:   Vec   *vecs;  /* holds the work vectors */
 33:   /* vv_allocated is the number of allocated gmres direction vectors */
 34:   PetscInt    q_preallocate,delta_allocate;
 35:   PetscInt    vv_allocated;
 36:   /* vecs_allocated is the total number of vecs available (used to 
 37:        simplify the dynamic allocation of vectors */
 38:   PetscInt    vecs_allocated;
 39:   /* Since we may call the user "obtain_work_vectors" several times, 
 40:        we have to keep track of the pointers that it has returned 
 41:       (so that we may free the storage) */
 42:   Vec         **user_work;
 43:   PetscInt    *mwork_alloc;    /* Number of work vectors allocated as part of
 44:                                a work-vector chunck */
 45:   PetscInt    nwork_alloc;     /* Number of work vectors allocated */

 47:   /* In order to allow the solution to be constructed during the solution
 48:      process, we need some additional information: */

 50:   PetscInt    it;              /* Current itethreshn: inside restart */
 51:   PetscScalar *nrs;            /* temp that holds the coefficients of the 
 52:                                Krylov vectors that form the minimum residual
 53:                                solution */
 54:   Vec         sol_temp;        /* used to hold temporary solution */
 55: 
 56:   /* Data specific to DGMRES */
 57:   Vec                        *U;        /* Vectors that form the basis of the invariant subspace */
 58:   PetscScalar        *T;        /* T=U^T*M^{-1}*A*U */
 59:   PetscScalar        *TF;        /* The factors L and U from T = P*L*U */
 60:   PetscBLASInt                 *InvP;        /* Permutation Vector from the LU factorization of T */
 61:   PetscInt                neig;        /* number of eigenvalues to extract at each restart */
 62:   PetscInt                r;                /* current number of deflated eigenvalues */
 63:   PetscInt                 max_neig;        /* Maximum number of eigenvalues to deflate */
 64:   PetscReal         lambdaN;        /* modulus of the largest eigenvalue of A */
 65:   PetscReal                smv;         /* smaller multiple of the remaining allowed number of steps -- used for the adaptive strategy */
 66:   PetscInt                force; /* Force the use of the deflation at the restart */
 67:   PetscInt                matvecs; /* Total number of matrix-vectors */
 68:   PetscInt                GreatestEig; /* Extract the greatest eigenvalues instead */
 69: 
 70:   /* Work spaces */
 71:   Vec                        *mu;        /* Save the product M^{-1}AU */
 72:   PetscScalar        *Sr;         /* Schur vectors to extract */
 73:   Vec                        *X;         /* Schurs Vectors Sr projected to the entire space */
 74:   Vec                        *mx;         /* store the product M^{-1}*A*X */
 75:   PetscScalar        *umx;         /* store the product U^T*M^{-1}*A*X */
 76:   PetscScalar        *xmu;         /* store the product X^T*M^{-1}*A*U */
 77:   PetscScalar        *xmx;        /* store the product X^T*M^{-1}*A*X */
 78:   PetscScalar        *x1;         /* store the product U^T*x */
 79:   PetscScalar        *x2;         /* store the product U^T*x */
 80:   PetscScalar         *Sr2;         /* Schur vectors at the improvement step */
 81:   PetscScalar        *auau;         /* product of (M*A*U)^T*M*A*U */
 82:   PetscScalar        *auu;         /* product of (M*A*U)^T*U */
 83: 
 84:   PetscScalar         *work;         /* work space for LAPACK functions */
 85:   PetscInt                *iwork;        /* work space for LAPACK functions */
 86:   PetscReal                *orth;         /* Coefficients for the orthogonalization */
 87: 
 88:   PetscInt                improve; /* 0 = do not improve the eigenvalues; This is an experimental option */
 89: 
 90: } KSP_DGMRES;

 92: PetscLogEvent KSP_DGMRESComputeDeflationData, KSP_DGMRESApplyDeflation;
 93: #define HH(a,b)  (dgmres->hh_origin + (b)*(dgmres->max_k+2)+(a))
 94: #define HES(a,b) (dgmres->hes_origin + (b)*(dgmres->max_k+1)+(a))
 95: #define CC(a)    (dgmres->cc_origin + (a))
 96: #define SS(a)    (dgmres->ss_origin + (a))
 97: #define GRS(a)   (dgmres->rs_origin + (a))

 99: /* vector names */
100: #define VEC_OFFSET     2
101: #define VEC_TEMP       dgmres->vecs[0]
102: #define VEC_TEMP_MATOP dgmres->vecs[1]
103: #define VEC_VV(i)      dgmres->vecs[VEC_OFFSET+i]

105: #define EIG_OFFSET                        2
106: #define DGMRES_DEFAULT_EIG        1
107: #define DGMRES_DEFAULT_MAXEIG 100

109: #define        UU                dgmres->U
110: #define        TT                dgmres->T
111: #define        TTF                dgmres->TF
112: #define        XX                dgmres->X
113: #define        INVP        dgmres->InvP
114: #define        MU                dgmres->mu
115: #define        MX                dgmres->mx
116: #define        UMX                dgmres->umx
117: #define        XMU                dgmres->xmu
118: #define        XMX                dgmres->xmx
119: #define        X1                dgmres->x1
120: #define        X2                dgmres->x2
121: #define        SR                dgmres->Sr
122: #define        SR2                dgmres->Sr2
123: #define AUAU        dgmres->auau
124: #define AUU                dgmres->auu
125: #define        MAX_K        dgmres->max_k
126: #define        MAX_NEIG dgmres->max_neig
127: #define WORK        dgmres->work
128: #define        IWORK        dgmres->iwork
129: #define        ORTH        dgmres->orth
130: #define SMV 1
131: #endif