1: /*
2: Private data structure used by the GMRES method. This data structure
3: must be identical to the beginning of the KSP_FGMRES data structure
4: so if you CHANGE anything here you must also change it there.
5: */
9: #include <petsc/private/kspimpl.h> 11: #define KSPGMRESHEADER \ 12: /* Hessenberg matrix and orthogonalization information. */ \ 13: PetscScalar *hh_origin; /* holds hessenburg matrix that has been multiplied by plane rotations (upper tri) */ \ 14: PetscScalar *hes_origin; /* holds the original (unmodified) hessenberg matrix which may be used to estimate the Singular Values of the matrix */ \ 15: PetscScalar *hes_ritz; /* holds the last full Hessenberg matrix to compute (harmonic) Ritz pairs */ \ 16: PetscScalar *cc_origin; /* holds cosines for rotation matrices */ \ 17: PetscScalar *ss_origin; /* holds sines for rotation matrices */ \ 18: PetscScalar *rs_origin; /* holds the right-hand-side of the Hessenberg system */ \ 19: \ 20: PetscScalar *orthogwork; /* holds dot products computed in orthogonalization */ \ 21: \ 22: /* Work space for computing eigenvalues/singular values */ \ 23: PetscReal *Dsvd; \ 24: PetscScalar *Rsvd; \ 25: \ 26: \ 27: PetscReal haptol; /* tolerance for happy ending */ \ 28: PetscInt max_k; /* number of vectors in Krylov space, restart size */ \ 29: PetscInt nextra_vecs; /* number of extra vecs needed, e.g. for a pipeline */ \ 30: \ 31: PetscErrorCode (*orthog)(KSP,PetscInt); \ 32: KSPGMRESCGSRefinementType cgstype; \ 33: \ 34: Vec *vecs; /* the work vectors */ \ 35: Vec *vecb; /* holds the last full basis vectors of the Krylov subspace to compute (harmonic) Ritz pairs */ \ 36: PetscInt q_preallocate; /* 0=don't preallocate space for work vectors */ \ 37: PetscInt delta_allocate; /* number of vectors to preallocaate in each block if not preallocated */ \ 38: PetscInt vv_allocated; /* number of allocated gmres direction vectors */ \ 39: PetscInt vecs_allocated; /* total number of vecs available */ \ 40: /* Since we may call the user "obtain_work_vectors" several times, we have to keep track of the pointers that it has returned */ \ 41: Vec **user_work; \ 42: PetscInt *mwork_alloc; /* Number of work vectors allocated as part of a work-vector chunck */ \ 43: PetscInt nwork_alloc; /* Number of work vector chunks allocated */ \ 44: \ 45: /* Information for building solution */ \ 46: PetscInt it; /* Current iteration: inside restart */ \ 47: PetscInt fullcycle; /* Current number of complete cycle */ \ 48: PetscScalar *nrs; /* temp that holds the coefficients of the Krylov vectors that form the minimum residual solution */ \ 49: Vec sol_temp; /* used to hold temporary solution */ 51: typedef struct {
52: KSPGMRESHEADER
53: } KSP_GMRES;
55: PETSC_INTERN PetscErrorCode KSPView_GMRES(KSP,PetscViewer);
56: PETSC_INTERN PetscErrorCode KSPSetUp_GMRES(KSP);
57: PETSC_INTERN PetscErrorCode KSPSetFromOptions_GMRES(PetscOptionItems *PetscOptionsObject,KSP);
58: PETSC_INTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal*,PetscReal*);
59: PETSC_INTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*);
60: PETSC_INTERN PetscErrorCode KSPComputeRitz_GMRES(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal*,PetscReal*);
61: PETSC_INTERN PetscErrorCode KSPReset_GMRES(KSP);
62: PETSC_INTERN PetscErrorCode KSPDestroy_GMRES(KSP);
63: PETSC_INTERN PetscErrorCode KSPGMRESGetNewVectors(KSP,PetscInt);
65: typedef PetscErrorCode (*FCN)(KSP,PetscInt); /* force argument to next function to not be extern C*/
67: PETSC_INTERN PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP,PetscReal);
68: PETSC_INTERN PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP);
69: PETSC_INTERN PetscErrorCode KSPGMRESSetRestart_GMRES(KSP,PetscInt);
70: PETSC_INTERN PetscErrorCode KSPGMRESGetRestart_GMRES(KSP,PetscInt*);
71: PETSC_INTERN PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP,FCN);
72: PETSC_INTERN PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP,FCN*);
73: PETSC_INTERN PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);
74: PETSC_INTERN PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType*);
76: /* These macros are guarded because they are redefined by derived implementations */
77: #if !defined(KSPGMRES_NO_MACROS)
78: #define HH(a,b) (gmres->hh_origin + (b)*(gmres->max_k+2)+(a)) 79: #define HES(a,b) (gmres->hes_origin + (b)*(gmres->max_k+1)+(a)) 80: #define CC(a) (gmres->cc_origin + (a)) 81: #define SS(a) (gmres->ss_origin + (a)) 82: #define GRS(a) (gmres->rs_origin + (a)) 84: /* vector names */
85: #define VEC_OFFSET 2 86: #define VEC_TEMP gmres->vecs[0] 87: #define VEC_TEMP_MATOP gmres->vecs[1] 88: #define VEC_VV(i) gmres->vecs[VEC_OFFSET+i] 89: #endif
91: #endif