Actual source code: gmresp.h

  1: /* $Id: gmresp.h,v 1.18 2001/08/07 03:03:51 balay Exp $ */
  2: /*
  3:    Private data structure used by the GMRES method.
  4: */

 8:  #include src/sles/ksp/kspimpl.h

 10: typedef struct {
 11:     /* Hessenberg matrix and orthogonalization information.  Hes holds
 12:        the original (unmodified) hessenberg matrix which may be used
 13:        to estimate the Singular Values of the matrix */
 14:     PetscScalar *hh_origin,*hes_origin,*cc_origin,*ss_origin,*rs_origin;
 15:     /* Work space for computing eigenvalues/singular values */
 16:     PetscReal *Dsvd;
 17:     PetscScalar *Rsvd;
 18: 
 19:     /* parameters */
 20:     PetscReal haptol;
 21:     int    max_k;

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

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

 43:     int    it;              /* Current iteration: inside restart */
 44:     PetscScalar *nrs;            /* temp that holds the coefficients of the 
 45:                                Krylov vectors that form the minimum residual
 46:                                solution */
 47:     Vec    sol_temp;        /* used to hold temporary solution */

 49: } KSP_GMRES;

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

 57: /* vector names */
 58: #define VEC_OFFSET     2
 59: #define VEC_SOLN       ksp->vec_sol
 60: #define VEC_RHS        ksp->vec_rhs
 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