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