Actual source code: lgmresp.h

  1: /* A. Baker */
  2: /*
  3:    Private data structure used by the LGMRES method.
  4: */


 9:  #include src/ksp/ksp/kspimpl.h

 11:   typedef struct {
 12:     /* Hessenberg matrix and orthogonalization information. */
 13:     PetscScalar *hh_origin;   /* holds hessenburg matrix that has been
 14:                             multiplied by plane rotations (upper tri) */
 15:     PetscScalar *hes_origin;  /* holds the original (unmodified) hessenberg matrix 
 16:                             which may be used to estimate the Singular Values
 17:                             of the matrix */
 18:     PetscScalar *cc_origin;   /* holds cosines for rotation matrices */
 19:     PetscScalar *ss_origin;   /* holds sines for rotation matrices */
 20:     PetscScalar *rs_origin;   /* holds the right-hand-side of the Hessenberg system */

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

 31:     PetscErrorCode (*orthog)(KSP,PetscInt); /* orthogonalization function to use */
 32:     KSPGMRESCGSRefinementType cgstype;
 33: 
 34:     Vec         *vecs;              /* holds the work vectors */
 35: 
 36:     PetscInt    q_preallocate;     /* 0 = don't pre-allocate space for work vectors */
 37:     PetscInt    delta_allocate;    /* the number of vectors to allocate in each block 
 38:                                       if not pre-allocated */
 39:     PetscInt    vv_allocated;      /* vv_allocated is the number of allocated lgmres 
 40:                                       direction vectors */
 41: 
 42:     PetscInt    vecs_allocated;    /* vecs_allocated is the total number of vecs 
 43:                                       available - used to simplify the dynamic
 44:                                      allocation of vectors */
 45: 
 46:     Vec         **user_work;       /* Since we may call the user "obtain_work_vectors" 
 47:                                      several times, we have to keep track of the pointers
 48:                                      that it has returned (so that we may free the 
 49:                                      storage) */

 51:     PetscInt    *mwork_alloc;      /* Number of work vectors allocated as part of
 52:                                       a work-vector chunck */
 53:     PetscInt    nwork_alloc;       /* Number of work-vector chunks allocated */


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

 59:     PetscInt    it;              /* Current iteration */
 60:     PetscScalar *nrs;            /* temp that holds the coefficients of the 
 61:                                     Krylov vectors that form the minimum residual
 62:                                     solution */
 63:     Vec         sol_temp;        /* used to hold temporary solution */


 66:     /* LGMRES_MOD - make these for the z vectors - new storage for lgmres */
 67:     Vec         *augvecs;            /* holds the error approximation vectors for lgmres. */
 68:     Vec         **augvecs_user_work; /* same purpose as user_work above, but this one is
 69:                                          for our error approx vectors */
 70:     PetscInt    aug_vv_allocated;      /* aug_vv_allocated is the number of allocated lgmres 
 71:                                           augmentation vectors */
 72:     PetscInt    aug_vecs_allocated;    /* aug_vecs_allocated is the total number of augmentation vecs 
 73:                                           available - used to simplify the dynamic
 74:                                        allocation of vectors */

 76:     PetscInt    augwork_alloc;       /*size of chunk allocated for augmentation vectors */

 78:     PetscInt    aug_dim;             /* max number of augmented directions to add */

 80:     PetscInt    aug_ct;              /* number of aug. vectors available */

 82:     PetscInt    *aug_order;          /*keeps track of order to use aug. vectors*/

 84:     PetscInt    approx_constant;   /* = 1 then the approx space at each restart will 
 85:                                   be  size max_k .  Therefore, more than (max_k - aug_dim) 
 86:                                   krylov vectors may be used if less than aug_dim error 
 87:                                   approximations are available (in the first few restarts,
 88:                                   for example) to keep the space a constant size. */
 89: 
 90:     PetscInt    matvecs;            /*keep track of matvecs */
 91: } KSP_LGMRES;


 94: #define HH(a,b)  (lgmres->hh_origin + (b)*(lgmres->max_k+2)+(a)) 
 95:                  /* HH will be size (max_k+2)*(max_k+1)  -  think of HH as 
 96:                     being stored columnwise (inc. zeros) for access purposes. */
 97: #define HES(a,b) (lgmres->hes_origin + (b)*(lgmres->max_k+1)+(a)) 
 98:                   /* HES will be size (max_k + 1) * (max_k + 1) - 
 99:                      again, think of HES as being stored columnwise */
100: #define CC(a)    (lgmres->cc_origin + (a)) /* CC will be length (max_k+1) - cosines */
101: #define SS(a)    (lgmres->ss_origin + (a)) /* SS will be length (max_k+1) - sines */
102: #define GRS(a)    (lgmres->rs_origin + (a)) /* GRS will be length (max_k+2) - rt side */

104: /* vector names */
105: #define VEC_OFFSET     2
106: #define VEC_TEMP       lgmres->vecs[0]               /* work space */  
107: #define VEC_TEMP_MATOP lgmres->vecs[1]               /* work space */
108: #define VEC_VV(i)      lgmres->vecs[VEC_OFFSET+i]    /* use to access
109:                                                         othog basis vectors */
110: /*LGMRES_MOD */
111: #define AUG_OFFSET     1
112: #define AUGVEC(i)      lgmres->augvecs[AUG_OFFSET+i]   /*error approx vecors */
113: #define AUG_ORDER(i)   lgmres->aug_order[i]            /*order in which to augment */ 
114: #define A_AUGVEC(i)    lgmres->augvecs[AUG_OFFSET+i+lgmres->aug_dim] /*A times error vector */
115: #define AUG_TEMP       lgmres->augvecs[0]              /* work vector */ 
116: #endif