Actual source code: lgmresimpl.h
1: /* A. Baker */
2: /*
3: Private data structure used by the LGMRES method.
4: */
9: #define KSPGMRES_NO_MACROS
10: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
12: typedef struct {
13: KSPGMRESHEADER
15: /* LGMRES_MOD - make these for the z vectors - new storage for lgmres */
16: Vec *augvecs; /* holds the error approximation vectors for lgmres. */
17: Vec **augvecs_user_work; /* same purpose as user_work above, but this one is
18: for our error approx vectors */
19: /* currently only augvecs_user_work[0] is used, not sure if this will be */
20: /* extended in the future to use more, or if this is a design bug */
21: PetscInt aug_vv_allocated; /* aug_vv_allocated is the number of allocated lgmres
22: augmentation vectors */
23: PetscInt aug_vecs_allocated; /* aug_vecs_allocated is the total number of augmentation vecs
24: available - used to simplify the dynamic
25: allocation of vectors */
26: PetscScalar *hwork; /* work array to hold Hessenberg product */
28: PetscInt augwork_alloc; /*size of chunk allocated for augmentation vectors */
30: PetscInt aug_dim; /* max number of augmented directions to add */
32: PetscInt aug_ct; /* number of aug. vectors available */
34: PetscInt *aug_order; /*keeps track of order to use aug. vectors*/
36: PetscBool approx_constant; /* = 1 then the approx space at each restart will
37: be size max_k . Therefore, more than (max_k - aug_dim)
38: krylov vectors may be used if less than aug_dim error
39: approximations are available (in the first few restarts,
40: for example) to keep the space a constant size. */
42: PetscInt matvecs; /*keep track of matvecs */
43: } KSP_LGMRES;
45: #define HH(a,b) (lgmres->hh_origin + (b)*(lgmres->max_k+2)+(a))
46: /* HH will be size (max_k+2)*(max_k+1) - think of HH as
47: being stored columnwise (inc. zeros) for access purposes. */
48: #define HES(a,b) (lgmres->hes_origin + (b)*(lgmres->max_k+1)+(a))
49: /* HES will be size (max_k + 1) * (max_k + 1) -
50: again, think of HES as being stored columnwise */
51: #define CC(a) (lgmres->cc_origin + (a)) /* CC will be length (max_k+1) - cosines */
52: #define SS(a) (lgmres->ss_origin + (a)) /* SS will be length (max_k+1) - sines */
53: #define GRS(a) (lgmres->rs_origin + (a)) /* GRS will be length (max_k+2) - rt side */
55: /* vector names */
56: #define VEC_OFFSET 2
57: #define VEC_TEMP lgmres->vecs[0] /* work space */
58: #define VEC_TEMP_MATOP lgmres->vecs[1] /* work space */
59: #define VEC_VV(i) lgmres->vecs[VEC_OFFSET+i] /* use to access
60: othog basis vectors */
61: /*LGMRES_MOD */
62: #define AUG_OFFSET 1
63: #define AUGVEC(i) lgmres->augvecs[AUG_OFFSET+i] /*error approx vectors */
64: #define AUG_ORDER(i) lgmres->aug_order[i] /*order in which to augment */
65: #define A_AUGVEC(i) lgmres->augvecs[AUG_OFFSET+i+lgmres->aug_dim] /*A times error vector */
66: #define AUG_TEMP lgmres->augvecs[0] /* work vector */
67: #endif