Actual source code: lgmresimpl.h

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

  5: #pragma once

  7: #define KSPGMRES_NO_MACROS
  8: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>

 10: typedef struct {
 11:   KSPGMRESHEADER

 13:   /* LGMRES_MOD - make these for the z vectors - new storage for lgmres */
 14:   Vec  *augvecs;           /* holds the error approximation vectors for lgmres. */
 15:   Vec **augvecs_user_work; /* same purpose as user_work above, but this one is for our error approx vectors */
 16:   /* currently only augvecs_user_work[0] is used, not sure if this will be */
 17:   /* extended in the future to use more, or if this is a design bug */
 18:   PetscInt     aug_vv_allocated;   /* aug_vv_allocated is the number of allocated lgmres augmentation vectors */
 19:   PetscInt     aug_vecs_allocated; /* aug_vecs_allocated is the total number of augmentation vecs available */
 20:   PetscScalar *hwork;              /* work array to hold Hessenberg product */

 22:   PetscInt  augwork_alloc;   /*size of chunk allocated for augmentation vectors */
 23:   PetscInt  aug_dim;         /* max number of augmented directions to add */
 24:   PetscInt  aug_ct;          /* number of aug. vectors available */
 25:   PetscInt *aug_order;       /*keeps track of order to use aug. vectors*/
 26:   PetscBool approx_constant; /* = 1 then the approx space at each restart will
 27:                                   be  size max_k .  Therefore, more than (max_k - aug_dim)
 28:                                   krylov vectors may be used if less than aug_dim error
 29:                                   approximations are available (in the first few restarts,
 30:                                   for example) to keep the space a constant size. */

 32:   PetscInt matvecs; /*keep track of matvecs */
 33: } KSP_LGMRES;

 35: #define HH(a, b) (lgmres->hh_origin + (b) * (lgmres->max_k + 2) + (a))
 36: /* HH will be size (max_k+2)*(max_k+1)  -  think of HH as
 37:    being stored columnwise (inc. zeros) for access purposes. */
 38: #define HES(a, b) (lgmres->hes_origin + (b) * (lgmres->max_k + 1) + (a))
 39: /* HES will be size (max_k + 1) * (max_k + 1) -
 40:    again, think of HES as being stored columnwise */
 41: #define CC(a)  (lgmres->cc_origin + (a)) /* CC will be length (max_k+1) - cosines */
 42: #define SS(a)  (lgmres->ss_origin + (a)) /* SS will be length (max_k+1) - sines */
 43: #define GRS(a) (lgmres->rs_origin + (a)) /* GRS will be length (max_k+2) - rt side */

 45: /* vector names */
 46: #define VEC_OFFSET     2
 47: #define VEC_TEMP       lgmres->vecs[0]              /* work space */
 48: #define VEC_TEMP_MATOP lgmres->vecs[1]              /* work space */
 49: #define VEC_VV(i)      lgmres->vecs[VEC_OFFSET + i] /* use to access othog basis vectors */
 50: /*LGMRES_MOD */
 51: #define AUG_OFFSET   1
 52: #define AUGVEC(i)    lgmres->augvecs[AUG_OFFSET + i]                   /*error approx vectors */
 53: #define AUG_ORDER(i) lgmres->aug_order[i]                              /*order in which to augment */
 54: #define A_AUGVEC(i)  lgmres->augvecs[AUG_OFFSET + i + lgmres->aug_dim] /*A times error vector */
 55: #define AUG_TEMP     lgmres->augvecs[0]                                /* work vector */

 57: #define LGMRES_DELTA_DIRECTIONS 10
 58: #define LGMRES_DEFAULT_MAXK     30
 59: #define LGMRES_DEFAULT_AUGDIM   2 /*default number of augmentation vectors */