Actual source code: gmresp.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 src/ksp/ksp/kspimpl.h

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

 17:   /* Work space for computing eigenvalues/singular values */
 18:   PetscReal   *Dsvd;
 19:   PetscScalar *Rsvd;
 20: 
 21:   /* parameters */
 22:   PetscReal haptol;
 23:   PetscInt  max_k;

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

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

 46:   PetscInt    it;              /* Current iteration: inside restart */
 47:   PetscScalar *nrs;            /* temp that holds the coefficients of the 
 48:                                Krylov vectors that form the minimum residual
 49:                                solution */
 50:   Vec         sol_temp;        /* used to hold temporary solution */
 51: } KSP_GMRES;

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

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

 65: #endif