Actual source code: fgmresp.h

  1: /*
  2:    Private data structure used by the FGMRES method. The beginning of this
  3:  data structure MUST be identical to the beginning of KSP_GMRES since they
  4:  share several functions!
  5: */


 10:  #include include/private/kspimpl.h

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

 23:     PetscScalar *orthogwork; /* holds dot products computed in orthogonalization */

 25:     /* Work space for computing eigenvalues/singular values */
 26:     PetscReal   *Dsvd;
 27:     PetscScalar *Rsvd;
 28: 
 29:     /* parameters */
 30:     PetscReal haptol;            /* tolerance used for the "HAPPY BREAK DOWN"  */
 31:     PetscInt  max_k;             /* maximum number of Krylov directions to find 
 32:                                     before restarting */

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

 54:     PetscInt *mwork_alloc;       /* Number of work vectors allocated as part of
 55:                                     a work-vector chunck */
 56:     PetscInt nwork_alloc;         /* Number of work-vector chunks allocated */


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

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


 69:     /* new storage for fgmres */
 70:     Vec         *prevecs;        /* holds the preconditioned basis vectors for fgmres.  
 71:                                     We will allocate these at the same time as vecs 
 72:                                     above (and in the same "chunks". */
 73:     Vec          **prevecs_user_work; /* same purpose as user_work above, but this one is
 74:                                     for our preconditioned vectors */

 76:     /* we need a function for interacting with the pcfamily */
 77: 
 78:     PetscErrorCode (*modifypc)(KSP,PetscInt,PetscInt,PetscReal,void*);  /* function to modify the preconditioner*/
 79:     PetscErrorCode (*modifydestroy)(void*);
 80:     void   *modifyctx;
 81: } KSP_FGMRES;

 83: #define HH(a,b)  (fgmres->hh_origin + (b)*(fgmres->max_k+2)+(a)) 
 84:                  /* HH will be size (max_k+2)*(max_k+1)  -  think of HH as 
 85:                     being stored columnwise for access purposes. */
 86: #define HES(a,b) (fgmres->hes_origin + (b)*(fgmres->max_k+1)+(a)) 
 87:                   /* HES will be size (max_k + 1) * (max_k + 1) - 
 88:                      again, think of HES as being stored columnwise */
 89: #define CC(a)    (fgmres->cc_origin + (a)) /* CC will be length (max_k+1) - cosines */
 90: #define SS(a)    (fgmres->ss_origin + (a)) /* SS will be length (max_k+1) - sines */
 91: #define RS(a)    (fgmres->rs_origin + (a)) /* RS will be length (max_k+2) - rt side */

 93: /* vector names */
 94: #define VEC_OFFSET     2
 95: #define VEC_TEMP       fgmres->vecs[0]               /* work space */  
 96: #define VEC_TEMP_MATOP fgmres->vecs[1]               /* work space */
 97: #define VEC_VV(i)      fgmres->vecs[VEC_OFFSET+i]    /* use to access
 98:                                                         othog basis vectors */
 99: #define PREVEC(i)      fgmres->prevecs[i]            /* use to access 
100:                                                         preconditioned basis */

102: #endif